1807171c5SHong Zhang 2*3b1b9624SHong Zhang /* 3*3b1b9624SHong Zhang Defines projective product routines where A is a SeqAIJ matrix 4*3b1b9624SHong Zhang C = R * A * R^T 5*3b1b9624SHong Zhang */ 6*3b1b9624SHong 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 12*3b1b9624SHong 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); 26*3b1b9624SHong 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__ 38807171c5SHong Zhang #define __FUNCT__ "MatRARtSymbolic_SeqAIJ_SeqAIJ" 39807171c5SHong Zhang PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ(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; 48*3b1b9624SHong Zhang #if defined(RART_PROFILE) 49807171c5SHong Zhang PetscLogDouble GColor=0.0,MCCreate=0.0,MDenCreate=0.0,t0,tf,etime=0.0; 50*3b1b9624SHong Zhang #endif 51807171c5SHong Zhang Mat_SeqAIJ *c; 52807171c5SHong Zhang 53807171c5SHong Zhang PetscFunctionBegin; 54*3b1b9624SHong Zhang #if defined(RART_PROFILE) 558563dfccSBarry Smith ierr = PetscTime(&t0);CHKERRQ(ierr); 56*3b1b9624SHong 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; 65807171c5SHong Zhang 66807171c5SHong Zhang /* create a supporting struct */ 67807171c5SHong Zhang ierr = PetscNew(Mat_RARt,&rart);CHKERRQ(ierr); 68257c235dSHong Zhang c = (Mat_SeqAIJ*)(*C)->data; 69257c235dSHong Zhang c->rart = rart; 70*3b1b9624SHong Zhang #if defined(RART_PROFILE) 718563dfccSBarry Smith ierr = PetscTime(&tf);CHKERRQ(ierr); 72807171c5SHong Zhang etime += tf - t0; 73*3b1b9624SHong Zhang #endif 74807171c5SHong Zhang 7522e94b5dSHong Zhang /* ------ Use coloring ---------- */ 76807171c5SHong Zhang /* Create MatTransposeColoring from symbolic C=R*A*R^T */ 77*3b1b9624SHong Zhang #if defined(RART_PROFILE) 788563dfccSBarry Smith ierr = PetscTime(&t0);CHKERRQ(ierr); 79*3b1b9624SHong Zhang #endif 80807171c5SHong Zhang ierr = MatGetColoring(*C,MATCOLORINGLF,&iscoloring);CHKERRQ(ierr); 81*3b1b9624SHong Zhang #if defined(RART_PROFILE) 828563dfccSBarry Smith ierr = PetscTime(&tf);CHKERRQ(ierr); 83807171c5SHong Zhang GColor += tf - t0; 848563dfccSBarry Smith ierr = PetscTime(&t0);CHKERRQ(ierr); 85*3b1b9624SHong Zhang #endif 86807171c5SHong Zhang ierr = MatTransposeColoringCreate(*C,iscoloring,&matcoloring);CHKERRQ(ierr); 872205254eSKarl Rupp 88807171c5SHong Zhang rart->matcoloring = matcoloring; 892205254eSKarl Rupp 90807171c5SHong Zhang ierr = ISColoringDestroy(&iscoloring);CHKERRQ(ierr); 91*3b1b9624SHong Zhang #if defined(RART_PROFILE) 928563dfccSBarry Smith ierr = PetscTime(&tf);CHKERRQ(ierr); 93807171c5SHong Zhang MCCreate += tf - t0; 948563dfccSBarry Smith ierr = PetscTime(&t0);CHKERRQ(ierr); 95*3b1b9624SHong Zhang #endif 96*3b1b9624SHong Zhang 97807171c5SHong Zhang /* Create Rt_dense */ 98807171c5SHong Zhang ierr = MatCreate(PETSC_COMM_SELF,&Rt_dense);CHKERRQ(ierr); 99807171c5SHong Zhang ierr = MatSetSizes(Rt_dense,A->cmap->n,matcoloring->ncolors,A->cmap->n,matcoloring->ncolors);CHKERRQ(ierr); 100807171c5SHong Zhang ierr = MatSetType(Rt_dense,MATSEQDENSE);CHKERRQ(ierr); 1010298fd71SBarry Smith ierr = MatSeqDenseSetPreallocation(Rt_dense,NULL);CHKERRQ(ierr); 1022205254eSKarl Rupp 103807171c5SHong Zhang Rt_dense->assembled = PETSC_TRUE; 104807171c5SHong Zhang rart->Rt = Rt_dense; 105807171c5SHong Zhang 106807171c5SHong Zhang /* Create RARt_dense = R*A*Rt_dense */ 107807171c5SHong Zhang ierr = MatCreate(PETSC_COMM_SELF,&RARt_dense);CHKERRQ(ierr); 108807171c5SHong Zhang ierr = MatSetSizes(RARt_dense,(*C)->rmap->n,matcoloring->ncolors,(*C)->rmap->n,matcoloring->ncolors);CHKERRQ(ierr); 109807171c5SHong Zhang ierr = MatSetType(RARt_dense,MATSEQDENSE);CHKERRQ(ierr); 1100298fd71SBarry Smith ierr = MatSeqDenseSetPreallocation(RARt_dense,NULL);CHKERRQ(ierr); 1112205254eSKarl Rupp 112807171c5SHong Zhang rart->RARt = RARt_dense; 113807171c5SHong Zhang 114807171c5SHong Zhang /* Allocate work array to store columns of A*R^T used in MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense() */ 115807171c5SHong Zhang ierr = PetscMalloc(A->rmap->n*4*sizeof(PetscScalar),&rart->work);CHKERRQ(ierr); 116807171c5SHong Zhang 117*3b1b9624SHong Zhang #if defined(RART_PROFILE) 1188563dfccSBarry Smith ierr = PetscTime(&tf);CHKERRQ(ierr); 119807171c5SHong Zhang MDenCreate += tf - t0; 120*3b1b9624SHong Zhang #endif 121807171c5SHong Zhang 122807171c5SHong Zhang rart->destroy = (*C)->ops->destroy; 123807171c5SHong Zhang (*C)->ops->destroy = MatDestroy_SeqAIJ_RARt; 124807171c5SHong Zhang 125807171c5SHong Zhang /* clean up */ 126807171c5SHong Zhang ierr = MatRestoreSymbolicTranspose_SeqAIJ(R,&rti,&rtj);CHKERRQ(ierr); 127807171c5SHong Zhang ierr = MatDestroy(&P);CHKERRQ(ierr); 128807171c5SHong Zhang 129807171c5SHong Zhang #if defined(PETSC_USE_INFO) 1301ce71dffSSatish Balay { 131807171c5SHong Zhang PetscReal density= (PetscReal)(c->nz)/(RARt_dense->rmap->n*RARt_dense->cmap->n); 132a2ea699eSBarry 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); 133*3b1b9624SHong Zhang #if defined(RART_PROFILE) 134a2ea699eSBarry 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); 135*3b1b9624SHong Zhang #endif 1361ce71dffSSatish Balay } 137807171c5SHong Zhang #endif 138807171c5SHong Zhang PetscFunctionReturn(0); 139807171c5SHong Zhang } 140807171c5SHong Zhang 141807171c5SHong Zhang /* 142807171c5SHong Zhang RAB = R * A * B, R and A in seqaij format, B in dense format; 143807171c5SHong Zhang */ 144807171c5SHong Zhang #undef __FUNCT__ 145807171c5SHong Zhang #define __FUNCT__ "MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense" 146807171c5SHong Zhang PetscErrorCode MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense(Mat R,Mat A,Mat B,Mat RAB,PetscScalar *work) 147807171c5SHong Zhang { 148807171c5SHong Zhang Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*r=(Mat_SeqAIJ*)R->data; 149807171c5SHong Zhang PetscErrorCode ierr; 150807171c5SHong Zhang PetscScalar *b,r1,r2,r3,r4,*b1,*b2,*b3,*b4; 151807171c5SHong Zhang MatScalar *aa,*ra; 152807171c5SHong Zhang PetscInt cn =B->cmap->n,bm=B->rmap->n,col,i,j,n,*ai=a->i,*aj,am=A->rmap->n; 153807171c5SHong Zhang PetscInt am2=2*am,am3=3*am,bm4=4*bm; 154807171c5SHong Zhang PetscScalar *d,*c,*c2,*c3,*c4; 155807171c5SHong Zhang PetscInt *rj,rm=R->rmap->n,dm=RAB->rmap->n,dn=RAB->cmap->n; 156807171c5SHong Zhang PetscInt rm2=2*rm,rm3=3*rm,colrm; 157*3b1b9624SHong Zhang #if defined(RART_PROFILE) 158274010c0SHong Zhang PetscLogDouble t0,tf,Mult_sp_den1=0.0,Mult_sp_den2=0.0; 159*3b1b9624SHong Zhang #endif 160807171c5SHong Zhang 161807171c5SHong Zhang PetscFunctionBegin; 162807171c5SHong Zhang if (!dm || !dn) PetscFunctionReturn(0); 163807171c5SHong 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); 164807171c5SHong 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); 165807171c5SHong 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); 166807171c5SHong 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); 167807171c5SHong Zhang 168274010c0SHong Zhang { /* 169274010c0SHong Zhang This approach is not as good as original ones (will be removed later), but it reveals that 170274010c0SHong Zhang AB_den=A*B takes almost all execution time in R*A*B for src/ksp/ksp/examples/tutorials/ex56.c 171274010c0SHong Zhang */ 172274010c0SHong Zhang PetscBool via_matmatmult=PETSC_FALSE; 173274010c0SHong Zhang ierr = PetscOptionsGetBool(NULL,"-matrart_via_matmatmult",&via_matmatmult,NULL);CHKERRQ(ierr); 174274010c0SHong Zhang if (via_matmatmult) { 175274010c0SHong Zhang Mat AB_den; 176274010c0SHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqDense(A,B,0.0,&AB_den);CHKERRQ(ierr); 177274010c0SHong Zhang 178*3b1b9624SHong Zhang #if defined(RART_PROFILE) 17926be7272SHong Zhang ierr = PetscTime(&t0);CHKERRQ(ierr); 180*3b1b9624SHong Zhang #endif 181274010c0SHong Zhang ierr = MatMatMultNumeric_SeqAIJ_SeqDense(A,B,AB_den);CHKERRQ(ierr); 182*3b1b9624SHong Zhang #if defined(RART_PROFILE) 18326be7272SHong Zhang ierr = PetscTime(&tf);CHKERRQ(ierr); 184274010c0SHong Zhang Mult_sp_den1 += tf - t0; 18526be7272SHong Zhang ierr = PetscTime(&t0);CHKERRQ(ierr); 186*3b1b9624SHong Zhang #endif 187274010c0SHong Zhang ierr = MatMatMultNumeric_SeqAIJ_SeqDense(R,AB_den,RAB);CHKERRQ(ierr); 188*3b1b9624SHong Zhang #if defined(RART_PROFILE) 18926be7272SHong Zhang ierr = PetscTime(&tf);CHKERRQ(ierr); 190274010c0SHong Zhang Mult_sp_den2 += tf - t0; 191274010c0SHong 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); 192*3b1b9624SHong Zhang #endif 193274010c0SHong Zhang ierr = MatDestroy(&AB_den);CHKERRQ(ierr); 194274010c0SHong Zhang PetscFunctionReturn(0); 195274010c0SHong Zhang } 196274010c0SHong Zhang } 197274010c0SHong Zhang 1988c778c55SBarry Smith ierr = MatDenseGetArray(B,&b);CHKERRQ(ierr); 1998c778c55SBarry Smith ierr = MatDenseGetArray(RAB,&d);CHKERRQ(ierr); 200807171c5SHong Zhang b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm; 201807171c5SHong Zhang c = work; c2 = c + am; c3 = c2 + am; c4 = c3 + am; 202807171c5SHong Zhang for (col=0; col<cn-4; col += 4) { /* over columns of C */ 203807171c5SHong Zhang for (i=0; i<am; i++) { /* over rows of A in those columns */ 204807171c5SHong Zhang r1 = r2 = r3 = r4 = 0.0; 205807171c5SHong Zhang n = ai[i+1] - ai[i]; 206807171c5SHong Zhang aj = a->j + ai[i]; 207807171c5SHong Zhang aa = a->a + ai[i]; 208807171c5SHong Zhang for (j=0; j<n; j++) { 209807171c5SHong Zhang r1 += (*aa)*b1[*aj]; 210807171c5SHong Zhang r2 += (*aa)*b2[*aj]; 211807171c5SHong Zhang r3 += (*aa)*b3[*aj]; 212807171c5SHong Zhang r4 += (*aa++)*b4[*aj++]; 213807171c5SHong Zhang } 214807171c5SHong Zhang c[i] = r1; 215807171c5SHong Zhang c[am + i] = r2; 216807171c5SHong Zhang c[am2 + i] = r3; 217807171c5SHong Zhang c[am3 + i] = r4; 218807171c5SHong Zhang } 219807171c5SHong Zhang b1 += bm4; 220807171c5SHong Zhang b2 += bm4; 221807171c5SHong Zhang b3 += bm4; 222807171c5SHong Zhang b4 += bm4; 223807171c5SHong Zhang 224807171c5SHong Zhang /* RAB[:,col] = R*C[:,col] */ 225807171c5SHong Zhang colrm = col*rm; 226807171c5SHong Zhang for (i=0; i<rm; i++) { /* over rows of R in those columns */ 227807171c5SHong Zhang r1 = r2 = r3 = r4 = 0.0; 228807171c5SHong Zhang n = r->i[i+1] - r->i[i]; 229807171c5SHong Zhang rj = r->j + r->i[i]; 230807171c5SHong Zhang ra = r->a + r->i[i]; 231807171c5SHong Zhang for (j=0; j<n; j++) { 232807171c5SHong Zhang r1 += (*ra)*c[*rj]; 233807171c5SHong Zhang r2 += (*ra)*c2[*rj]; 234807171c5SHong Zhang r3 += (*ra)*c3[*rj]; 235807171c5SHong Zhang r4 += (*ra++)*c4[*rj++]; 236807171c5SHong Zhang } 237807171c5SHong Zhang d[colrm + i] = r1; 238807171c5SHong Zhang d[colrm + rm + i] = r2; 239807171c5SHong Zhang d[colrm + rm2 + i] = r3; 240807171c5SHong Zhang d[colrm + rm3 + i] = r4; 241807171c5SHong Zhang } 242807171c5SHong Zhang } 243807171c5SHong Zhang for (; col<cn; col++) { /* over extra columns of C */ 244807171c5SHong Zhang for (i=0; i<am; i++) { /* over rows of A in those columns */ 245807171c5SHong Zhang r1 = 0.0; 246807171c5SHong Zhang n = a->i[i+1] - a->i[i]; 247807171c5SHong Zhang aj = a->j + a->i[i]; 248807171c5SHong Zhang aa = a->a + a->i[i]; 249807171c5SHong Zhang for (j=0; j<n; j++) { 250807171c5SHong Zhang r1 += (*aa++)*b1[*aj++]; 251807171c5SHong Zhang } 252807171c5SHong Zhang c[i] = r1; 253807171c5SHong Zhang } 254807171c5SHong Zhang b1 += bm; 255807171c5SHong Zhang 256807171c5SHong Zhang for (i=0; i<rm; i++) { /* over rows of R in those columns */ 257807171c5SHong Zhang r1 = 0.0; 258807171c5SHong Zhang n = r->i[i+1] - r->i[i]; 259807171c5SHong Zhang rj = r->j + r->i[i]; 260807171c5SHong Zhang ra = r->a + r->i[i]; 261807171c5SHong Zhang for (j=0; j<n; j++) { 262807171c5SHong Zhang r1 += (*ra++)*c[*rj++]; 263807171c5SHong Zhang } 264807171c5SHong Zhang d[col*rm + i] = r1; 265807171c5SHong Zhang } 266807171c5SHong Zhang } 267807171c5SHong Zhang ierr = PetscLogFlops(cn*2.0*(a->nz + r->nz));CHKERRQ(ierr); 268807171c5SHong Zhang 2698c778c55SBarry Smith ierr = MatDenseRestoreArray(B,&b);CHKERRQ(ierr); 2708c778c55SBarry Smith ierr = MatDenseRestoreArray(RAB,&d);CHKERRQ(ierr); 271807171c5SHong Zhang ierr = MatAssemblyBegin(RAB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 272807171c5SHong Zhang ierr = MatAssemblyEnd(RAB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 273807171c5SHong Zhang PetscFunctionReturn(0); 274807171c5SHong Zhang } 275807171c5SHong Zhang 276807171c5SHong Zhang #undef __FUNCT__ 277807171c5SHong Zhang #define __FUNCT__ "MatRARtNumeric_SeqAIJ_SeqAIJ" 278807171c5SHong Zhang PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ(Mat A,Mat R,Mat C) 279807171c5SHong Zhang { 280807171c5SHong Zhang PetscErrorCode ierr; 281257c235dSHong Zhang Mat_SeqAIJ *c = (Mat_SeqAIJ*)C->data; 282257c235dSHong Zhang Mat_RARt *rart=c->rart; 283807171c5SHong Zhang MatTransposeColoring matcoloring; 2848d0a38e4SHong Zhang Mat Rt,RARt; 285*3b1b9624SHong Zhang #if defined(RART_PROFILE) 286807171c5SHong Zhang PetscLogDouble Mult_sp_den=0.0,app1=0.0,app2=0.0,t0,tf; 287*3b1b9624SHong Zhang #endif 288807171c5SHong Zhang 289807171c5SHong Zhang PetscFunctionBegin; 290807171c5SHong Zhang /* Get dense Rt by Apply MatTransposeColoring to R */ 291807171c5SHong Zhang matcoloring = rart->matcoloring; 292807171c5SHong Zhang Rt = rart->Rt; 2932205254eSKarl Rupp 294*3b1b9624SHong Zhang #if defined(RART_PROFILE) 2958563dfccSBarry Smith ierr = PetscTime(&t0);CHKERRQ(ierr); 296*3b1b9624SHong Zhang #endif 297807171c5SHong Zhang ierr = MatTransColoringApplySpToDen(matcoloring,R,Rt);CHKERRQ(ierr); 298*3b1b9624SHong Zhang #if defined(RART_PROFILE) 2998563dfccSBarry Smith ierr = PetscTime(&tf);CHKERRQ(ierr); 300807171c5SHong Zhang app1 += tf - t0; 301*3b1b9624SHong Zhang #endif 302807171c5SHong Zhang 30350647e95SHong Zhang /* Get dense RARt = R*A*Rt -- dominates! */ 304*3b1b9624SHong Zhang #if defined(RART_PROFILE) 3058563dfccSBarry Smith ierr = PetscTime(&t0);CHKERRQ(ierr); 306*3b1b9624SHong Zhang #endif 307807171c5SHong Zhang RARt = rart->RARt; 308807171c5SHong Zhang ierr = MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense(R,A,Rt,RARt,rart->work);CHKERRQ(ierr); 309*3b1b9624SHong Zhang #if defined(RART_PROFILE) 3108563dfccSBarry Smith ierr = PetscTime(&tf);CHKERRQ(ierr); 311807171c5SHong Zhang Mult_sp_den += tf - t0; 312*3b1b9624SHong Zhang #endif 313807171c5SHong Zhang 314807171c5SHong Zhang /* Recover C from C_dense */ 315*3b1b9624SHong Zhang #if defined(RART_PROFILE) 3168563dfccSBarry Smith ierr = PetscTime(&t0);CHKERRQ(ierr); 317*3b1b9624SHong Zhang #endif 318807171c5SHong Zhang ierr = MatTransColoringApplyDenToSp(matcoloring,RARt,C);CHKERRQ(ierr); 319*3b1b9624SHong Zhang #if defined(RART_PROFILE) 3208563dfccSBarry Smith ierr = PetscTime(&tf);CHKERRQ(ierr); 321807171c5SHong Zhang app2 += tf - t0; 322*3b1b9624SHong Zhang #endif 323807171c5SHong Zhang 324*3b1b9624SHong Zhang #if defined(RART_PROFILE) 325*3b1b9624SHong 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); 326807171c5SHong Zhang #endif 327807171c5SHong Zhang PetscFunctionReturn(0); 328807171c5SHong Zhang } 329807171c5SHong Zhang 330807171c5SHong Zhang #undef __FUNCT__ 331807171c5SHong Zhang #define __FUNCT__ "MatRARt_SeqAIJ_SeqAIJ" 332807171c5SHong Zhang PetscErrorCode MatRARt_SeqAIJ_SeqAIJ(Mat A,Mat R,MatReuse scall,PetscReal fill,Mat *C) 333807171c5SHong Zhang { 334807171c5SHong Zhang PetscErrorCode ierr; 335b3b4fc5aSHong Zhang PetscBool usecoloring = PETSC_FALSE,color_art; 336807171c5SHong Zhang 337807171c5SHong Zhang PetscFunctionBegin; 33822e94b5dSHong Zhang ierr = PetscOptionsGetBool(NULL,"-matrart_color",&usecoloring,NULL);CHKERRQ(ierr); 33922e94b5dSHong Zhang 3404fa85fdeSHong Zhang if (!usecoloring) { /* Rt=R^T, C=R*A*Rt - avoid inefficient sparse inner products */ 34122e94b5dSHong Zhang Mat Rt; 3424fa85fdeSHong Zhang Mat_SeqAIJ *c; 3434fa85fdeSHong Zhang Mat_RARt *rart; 3444fa85fdeSHong Zhang 3454fa85fdeSHong Zhang if (scall == MAT_INITIAL_MATRIX) { 3464fa85fdeSHong Zhang ierr = MatTranspose(R,MAT_INITIAL_MATRIX,&Rt);CHKERRQ(ierr); 34722e94b5dSHong Zhang ierr = MatMatMatMult(R,A,Rt,scall,fill,C);CHKERRQ(ierr); 3484fa85fdeSHong Zhang 3494fa85fdeSHong Zhang ierr = PetscNew(Mat_RARt,&rart);CHKERRQ(ierr); 3504fa85fdeSHong Zhang rart->Rt = Rt; 3514fa85fdeSHong Zhang c = (Mat_SeqAIJ*)(*C)->data; 3524fa85fdeSHong Zhang c->rart = rart; 3534fa85fdeSHong Zhang rart->destroy = (*C)->ops->destroy; 3544fa85fdeSHong Zhang (*C)->ops->destroy = MatDestroy_SeqAIJ_RARt; 3554fa85fdeSHong Zhang } else { 3564fa85fdeSHong Zhang c = (Mat_SeqAIJ*)(*C)->data; 3574fa85fdeSHong Zhang rart = c->rart; 3584fa85fdeSHong Zhang Rt = rart->Rt; 3594fa85fdeSHong Zhang ierr = MatTranspose(R,MAT_REUSE_MATRIX,&Rt);CHKERRQ(ierr); 3604fa85fdeSHong Zhang ierr = MatMatMatMult(R,A,Rt,scall,fill,C);CHKERRQ(ierr); 3614fa85fdeSHong Zhang } 36222e94b5dSHong Zhang PetscFunctionReturn(0); 36322e94b5dSHong Zhang } 36422e94b5dSHong Zhang 36522e94b5dSHong Zhang /* use coloring */ 366b3b4fc5aSHong Zhang /*--------------*/ 367b3b4fc5aSHong Zhang ierr = PetscOptionsGetBool(NULL,"-matrart_color_art",&color_art,NULL);CHKERRQ(ierr); 368b3b4fc5aSHong Zhang if (color_art) { /* apply coloring to A*R^T, not C = R*A*R^T */ 36936104f73SHong Zhang Mat ARt,RARt; 370*3b1b9624SHong Zhang #if defined(RART_PROFILE) 371b3b4fc5aSHong Zhang PetscLogDouble t0,t1,t2; 372*3b1b9624SHong Zhang #endif 373*3b1b9624SHong Zhang Mat_SeqAIJ *c; 374*3b1b9624SHong Zhang Mat_RARt *rart; 375*3b1b9624SHong Zhang 37636104f73SHong Zhang if (scall == MAT_INITIAL_MATRIX) { 37736104f73SHong Zhang ierr = PetscLogEventBegin(MAT_RARtSymbolic,A,R,0,0);CHKERRQ(ierr); 378b3b4fc5aSHong Zhang /* must use '-mat_no_inode' with '-matmattransmult_color 1' - do not knwo why? */ 37936104f73SHong Zhang ierr = MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(A,R,fill,&ARt);CHKERRQ(ierr); 38036104f73SHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(R,ARt,fill,&RARt);CHKERRQ(ierr); 381*3b1b9624SHong Zhang 382*3b1b9624SHong Zhang *C = RARt; 383*3b1b9624SHong Zhang ierr = PetscNew(Mat_RARt,&rart);CHKERRQ(ierr); 384*3b1b9624SHong Zhang c = (Mat_SeqAIJ*)(*C)->data; 385*3b1b9624SHong Zhang c->rart = rart; 386*3b1b9624SHong Zhang rart->ARt = ARt; 387*3b1b9624SHong Zhang rart->destroy = RARt->ops->destroy; 388*3b1b9624SHong Zhang RARt->ops->destroy = MatDestroy_SeqAIJ_RARt; 38936104f73SHong Zhang ierr = PetscLogEventEnd(MAT_RARtSymbolic,A,R,0,0);CHKERRQ(ierr); 39036104f73SHong Zhang } 39136104f73SHong Zhang ierr = PetscLogEventBegin(MAT_RARtNumeric,A,R,0,0);CHKERRQ(ierr); 392*3b1b9624SHong Zhang c = (Mat_SeqAIJ*)(*C)->data; 393*3b1b9624SHong Zhang rart = c->rart; 394*3b1b9624SHong Zhang ARt = rart->ARt; 395*3b1b9624SHong Zhang #if defined(RART_PROFILE) 39626be7272SHong Zhang ierr = PetscTime(&t0);CHKERRQ(ierr); 397*3b1b9624SHong Zhang #endif 398b3b4fc5aSHong Zhang ierr = MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(A,R,ARt);CHKERRQ(ierr); /* dominate! */ 399*3b1b9624SHong Zhang #if defined(RART_PROFILE) 40026be7272SHong Zhang ierr = PetscTime(&t1);CHKERRQ(ierr); 401*3b1b9624SHong Zhang #endif 402*3b1b9624SHong Zhang ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ(R,ARt,*C);CHKERRQ(ierr); 403*3b1b9624SHong Zhang #if defined(RART_PROFILE) 40426be7272SHong Zhang ierr = PetscTime(&t2);CHKERRQ(ierr); 405b3b4fc5aSHong Zhang printf(" matrart_color_art_num = %g + %g = %g\n",t1-t0,t2-t1,t2-t0); 406*3b1b9624SHong Zhang #endif 40736104f73SHong Zhang ierr = PetscLogEventEnd(MAT_RARtNumeric,A,R,0,0);CHKERRQ(ierr); 40836104f73SHong Zhang PetscFunctionReturn(0); 40936104f73SHong Zhang } 410b3b4fc5aSHong Zhang 411b3b4fc5aSHong Zhang if (scall == MAT_INITIAL_MATRIX) { /* apply coloring to C=R*A*R^T */ 41250647e95SHong Zhang ierr = PetscLogEventBegin(MAT_RARtSymbolic,A,R,0,0);CHKERRQ(ierr); 413807171c5SHong Zhang ierr = MatRARtSymbolic_SeqAIJ_SeqAIJ(A,R,fill,C);CHKERRQ(ierr); 41450647e95SHong Zhang ierr = PetscLogEventEnd(MAT_RARtSymbolic,A,R,0,0);CHKERRQ(ierr); 415807171c5SHong Zhang } 41650647e95SHong Zhang ierr = PetscLogEventBegin(MAT_RARtNumeric,A,R,0,0);CHKERRQ(ierr); 417807171c5SHong Zhang ierr = MatRARtNumeric_SeqAIJ_SeqAIJ(A,R,*C);CHKERRQ(ierr); 41850647e95SHong Zhang ierr = PetscLogEventEnd(MAT_RARtNumeric,A,R,0,0);CHKERRQ(ierr); 419807171c5SHong Zhang PetscFunctionReturn(0); 420807171c5SHong Zhang } 421