1807171c5SHong Zhang 2807171c5SHong Zhang #include <../src/mat/impls/aij/seq/aij.h> 3807171c5SHong Zhang #include <../src/mat/utils/freespace.h> 4807171c5SHong Zhang #include <../src/mat/impls/dense/seq/dense.h> /*I "petscmat.h" I*/ 58563dfccSBarry Smith #include <petsctime.h> 6807171c5SHong Zhang 7807171c5SHong Zhang /* 8807171c5SHong Zhang Defines projective product routines where A is a SeqAIJ matrix 9807171c5SHong Zhang C = R * A * R^T 10807171c5SHong Zhang */ 11807171c5SHong Zhang 12807171c5SHong Zhang #undef __FUNCT__ 13257c235dSHong Zhang #define __FUNCT__ "MatDestroy_SeqAIJ_RARt" 14257c235dSHong Zhang PetscErrorCode MatDestroy_SeqAIJ_RARt(Mat A) 15807171c5SHong Zhang { 16807171c5SHong Zhang PetscErrorCode ierr; 17257c235dSHong Zhang Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 18257c235dSHong Zhang Mat_RARt *rart = a->rart; 19807171c5SHong Zhang 20807171c5SHong Zhang PetscFunctionBegin; 21807171c5SHong Zhang ierr = MatTransposeColoringDestroy(&rart->matcoloring);CHKERRQ(ierr); 22807171c5SHong Zhang ierr = MatDestroy(&rart->Rt);CHKERRQ(ierr); 23807171c5SHong Zhang ierr = MatDestroy(&rart->RARt);CHKERRQ(ierr); 24807171c5SHong Zhang ierr = PetscFree(rart->work);CHKERRQ(ierr); 25807171c5SHong Zhang 26807171c5SHong Zhang A->ops->destroy = rart->destroy; 27807171c5SHong Zhang if (A->ops->destroy) { 28807171c5SHong Zhang ierr = (*A->ops->destroy)(A);CHKERRQ(ierr); 29807171c5SHong Zhang } 30257c235dSHong Zhang ierr = PetscFree(rart);CHKERRQ(ierr); 31807171c5SHong Zhang PetscFunctionReturn(0); 32807171c5SHong Zhang } 33807171c5SHong Zhang 34807171c5SHong Zhang #undef __FUNCT__ 35807171c5SHong Zhang #define __FUNCT__ "MatRARtSymbolic_SeqAIJ_SeqAIJ" 36807171c5SHong Zhang PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat R,PetscReal fill,Mat *C) 37807171c5SHong Zhang { 38807171c5SHong Zhang PetscErrorCode ierr; 39807171c5SHong Zhang Mat P; 40807171c5SHong Zhang PetscInt *rti,*rtj; 41807171c5SHong Zhang Mat_RARt *rart; 42807171c5SHong Zhang MatTransposeColoring matcoloring; 43807171c5SHong Zhang ISColoring iscoloring; 448d0a38e4SHong Zhang Mat Rt_dense,RARt_dense; 45807171c5SHong Zhang PetscLogDouble GColor=0.0,MCCreate=0.0,MDenCreate=0.0,t0,tf,etime=0.0; 46807171c5SHong Zhang Mat_SeqAIJ *c; 47807171c5SHong Zhang 48807171c5SHong Zhang PetscFunctionBegin; 498563dfccSBarry Smith ierr = PetscTime(&t0);CHKERRQ(ierr); 50807171c5SHong Zhang /* create symbolic P=Rt */ 51807171c5SHong Zhang ierr = MatGetSymbolicTranspose_SeqAIJ(R,&rti,&rtj);CHKERRQ(ierr); 520298fd71SBarry Smith ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,R->cmap->n,R->rmap->n,rti,rtj,NULL,&P);CHKERRQ(ierr); 53807171c5SHong Zhang 54807171c5SHong Zhang /* get symbolic C=Pt*A*P */ 55807171c5SHong Zhang ierr = MatPtAPSymbolic_SeqAIJ_SeqAIJ(A,P,fill,C);CHKERRQ(ierr); 56fd8e3dafSHong Zhang (*C)->rmap->bs = R->rmap->bs; 57fd8e3dafSHong Zhang (*C)->cmap->bs = R->rmap->bs; 58807171c5SHong Zhang 59807171c5SHong Zhang /* create a supporting struct */ 60807171c5SHong Zhang ierr = PetscNew(Mat_RARt,&rart);CHKERRQ(ierr); 61257c235dSHong Zhang c = (Mat_SeqAIJ*)(*C)->data; 62257c235dSHong Zhang c->rart = rart; 638563dfccSBarry Smith ierr = PetscTime(&tf);CHKERRQ(ierr); 64807171c5SHong Zhang etime += tf - t0; 65807171c5SHong Zhang 6622e94b5dSHong Zhang /* ------ Use coloring ---------- */ 67807171c5SHong Zhang /* Create MatTransposeColoring from symbolic C=R*A*R^T */ 688563dfccSBarry Smith ierr = PetscTime(&t0);CHKERRQ(ierr); 69807171c5SHong Zhang ierr = MatGetColoring(*C,MATCOLORINGLF,&iscoloring);CHKERRQ(ierr); 708563dfccSBarry Smith ierr = PetscTime(&tf);CHKERRQ(ierr); 71807171c5SHong Zhang GColor += tf - t0; 72807171c5SHong Zhang 738563dfccSBarry Smith ierr = PetscTime(&t0);CHKERRQ(ierr); 74807171c5SHong Zhang ierr = MatTransposeColoringCreate(*C,iscoloring,&matcoloring);CHKERRQ(ierr); 752205254eSKarl Rupp 76807171c5SHong Zhang rart->matcoloring = matcoloring; 772205254eSKarl Rupp 78807171c5SHong Zhang ierr = ISColoringDestroy(&iscoloring);CHKERRQ(ierr); 798563dfccSBarry Smith ierr = PetscTime(&tf);CHKERRQ(ierr); 80807171c5SHong Zhang MCCreate += tf - t0; 81807171c5SHong Zhang 828563dfccSBarry Smith ierr = PetscTime(&t0);CHKERRQ(ierr); 83807171c5SHong Zhang /* Create Rt_dense */ 84807171c5SHong Zhang ierr = MatCreate(PETSC_COMM_SELF,&Rt_dense);CHKERRQ(ierr); 85807171c5SHong Zhang ierr = MatSetSizes(Rt_dense,A->cmap->n,matcoloring->ncolors,A->cmap->n,matcoloring->ncolors);CHKERRQ(ierr); 86807171c5SHong Zhang ierr = MatSetType(Rt_dense,MATSEQDENSE);CHKERRQ(ierr); 870298fd71SBarry Smith ierr = MatSeqDenseSetPreallocation(Rt_dense,NULL);CHKERRQ(ierr); 882205254eSKarl Rupp 89807171c5SHong Zhang Rt_dense->assembled = PETSC_TRUE; 90807171c5SHong Zhang rart->Rt = Rt_dense; 91807171c5SHong Zhang 92807171c5SHong Zhang /* Create RARt_dense = R*A*Rt_dense */ 93807171c5SHong Zhang ierr = MatCreate(PETSC_COMM_SELF,&RARt_dense);CHKERRQ(ierr); 94807171c5SHong Zhang ierr = MatSetSizes(RARt_dense,(*C)->rmap->n,matcoloring->ncolors,(*C)->rmap->n,matcoloring->ncolors);CHKERRQ(ierr); 95807171c5SHong Zhang ierr = MatSetType(RARt_dense,MATSEQDENSE);CHKERRQ(ierr); 960298fd71SBarry Smith ierr = MatSeqDenseSetPreallocation(RARt_dense,NULL);CHKERRQ(ierr); 972205254eSKarl Rupp 98807171c5SHong Zhang rart->RARt = RARt_dense; 99807171c5SHong Zhang 100807171c5SHong Zhang /* Allocate work array to store columns of A*R^T used in MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense() */ 101807171c5SHong Zhang ierr = PetscMalloc(A->rmap->n*4*sizeof(PetscScalar),&rart->work);CHKERRQ(ierr); 102807171c5SHong Zhang 1038563dfccSBarry Smith ierr = PetscTime(&tf);CHKERRQ(ierr); 104807171c5SHong Zhang MDenCreate += tf - t0; 105807171c5SHong Zhang 106807171c5SHong Zhang rart->destroy = (*C)->ops->destroy; 107807171c5SHong Zhang (*C)->ops->destroy = MatDestroy_SeqAIJ_RARt; 108807171c5SHong Zhang 109807171c5SHong Zhang /* clean up */ 110807171c5SHong Zhang ierr = MatRestoreSymbolicTranspose_SeqAIJ(R,&rti,&rtj);CHKERRQ(ierr); 111807171c5SHong Zhang ierr = MatDestroy(&P);CHKERRQ(ierr); 112807171c5SHong Zhang 113807171c5SHong Zhang #if defined(PETSC_USE_INFO) 1141ce71dffSSatish Balay { 115807171c5SHong Zhang PetscReal density= (PetscReal)(c->nz)/(RARt_dense->rmap->n*RARt_dense->cmap->n); 116a2ea699eSBarry 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); 117a2ea699eSBarry 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); 1181ce71dffSSatish Balay } 119807171c5SHong Zhang #endif 120807171c5SHong Zhang PetscFunctionReturn(0); 121807171c5SHong Zhang } 122807171c5SHong Zhang 123807171c5SHong Zhang /* 124807171c5SHong Zhang RAB = R * A * B, R and A in seqaij format, B in dense format; 125807171c5SHong Zhang */ 126807171c5SHong Zhang #undef __FUNCT__ 127807171c5SHong Zhang #define __FUNCT__ "MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense" 128807171c5SHong Zhang PetscErrorCode MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense(Mat R,Mat A,Mat B,Mat RAB,PetscScalar *work) 129807171c5SHong Zhang { 130807171c5SHong Zhang Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*r=(Mat_SeqAIJ*)R->data; 131807171c5SHong Zhang PetscErrorCode ierr; 132807171c5SHong Zhang PetscScalar *b,r1,r2,r3,r4,*b1,*b2,*b3,*b4; 133807171c5SHong Zhang MatScalar *aa,*ra; 134807171c5SHong Zhang PetscInt cn =B->cmap->n,bm=B->rmap->n,col,i,j,n,*ai=a->i,*aj,am=A->rmap->n; 135807171c5SHong Zhang PetscInt am2=2*am,am3=3*am,bm4=4*bm; 136807171c5SHong Zhang PetscScalar *d,*c,*c2,*c3,*c4; 137807171c5SHong Zhang PetscInt *rj,rm=R->rmap->n,dm=RAB->rmap->n,dn=RAB->cmap->n; 138807171c5SHong Zhang PetscInt rm2=2*rm,rm3=3*rm,colrm; 139274010c0SHong Zhang PetscLogDouble t0,tf,Mult_sp_den1=0.0,Mult_sp_den2=0.0; 140807171c5SHong Zhang 141807171c5SHong Zhang PetscFunctionBegin; 142807171c5SHong Zhang if (!dm || !dn) PetscFunctionReturn(0); 143807171c5SHong 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); 144807171c5SHong 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); 145807171c5SHong 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); 146807171c5SHong 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); 147807171c5SHong Zhang 148274010c0SHong Zhang { /* 149274010c0SHong Zhang This approach is not as good as original ones (will be removed later), but it reveals that 150274010c0SHong Zhang AB_den=A*B takes almost all execution time in R*A*B for src/ksp/ksp/examples/tutorials/ex56.c 151274010c0SHong Zhang */ 152274010c0SHong Zhang PetscBool via_matmatmult=PETSC_FALSE; 153274010c0SHong Zhang ierr = PetscOptionsGetBool(NULL,"-matrart_via_matmatmult",&via_matmatmult,NULL);CHKERRQ(ierr); 154274010c0SHong Zhang if (via_matmatmult) { 155274010c0SHong Zhang Mat AB_den; 156274010c0SHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqDense(A,B,0.0,&AB_den);CHKERRQ(ierr); 157274010c0SHong Zhang 158274010c0SHong Zhang ierr = PetscGetTime(&t0);CHKERRQ(ierr); 159274010c0SHong Zhang ierr = MatMatMultNumeric_SeqAIJ_SeqDense(A,B,AB_den);CHKERRQ(ierr); 160274010c0SHong Zhang ierr = PetscGetTime(&tf);CHKERRQ(ierr); 161274010c0SHong Zhang Mult_sp_den1 += tf - t0; 162274010c0SHong Zhang ierr = PetscGetTime(&t0);CHKERRQ(ierr); 163274010c0SHong Zhang ierr = MatMatMultNumeric_SeqAIJ_SeqDense(R,AB_den,RAB);CHKERRQ(ierr); 164274010c0SHong Zhang ierr = PetscGetTime(&tf);CHKERRQ(ierr); 165274010c0SHong Zhang Mult_sp_den2 += tf - t0; 166274010c0SHong 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); 167274010c0SHong Zhang ierr = MatDestroy(&AB_den);CHKERRQ(ierr); 168274010c0SHong Zhang PetscFunctionReturn(0); 169274010c0SHong Zhang } 170274010c0SHong Zhang } 171274010c0SHong Zhang 1728c778c55SBarry Smith ierr = MatDenseGetArray(B,&b);CHKERRQ(ierr); 1738c778c55SBarry Smith ierr = MatDenseGetArray(RAB,&d);CHKERRQ(ierr); 174807171c5SHong Zhang b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm; 175807171c5SHong Zhang c = work; c2 = c + am; c3 = c2 + am; c4 = c3 + am; 176807171c5SHong Zhang for (col=0; col<cn-4; col += 4) { /* over columns of C */ 177807171c5SHong Zhang for (i=0; i<am; i++) { /* over rows of A in those columns */ 178807171c5SHong Zhang r1 = r2 = r3 = r4 = 0.0; 179807171c5SHong Zhang n = ai[i+1] - ai[i]; 180807171c5SHong Zhang aj = a->j + ai[i]; 181807171c5SHong Zhang aa = a->a + ai[i]; 182807171c5SHong Zhang for (j=0; j<n; j++) { 183807171c5SHong Zhang r1 += (*aa)*b1[*aj]; 184807171c5SHong Zhang r2 += (*aa)*b2[*aj]; 185807171c5SHong Zhang r3 += (*aa)*b3[*aj]; 186807171c5SHong Zhang r4 += (*aa++)*b4[*aj++]; 187807171c5SHong Zhang } 188807171c5SHong Zhang c[i] = r1; 189807171c5SHong Zhang c[am + i] = r2; 190807171c5SHong Zhang c[am2 + i] = r3; 191807171c5SHong Zhang c[am3 + i] = r4; 192807171c5SHong Zhang } 193807171c5SHong Zhang b1 += bm4; 194807171c5SHong Zhang b2 += bm4; 195807171c5SHong Zhang b3 += bm4; 196807171c5SHong Zhang b4 += bm4; 197807171c5SHong Zhang 198807171c5SHong Zhang /* RAB[:,col] = R*C[:,col] */ 199807171c5SHong Zhang colrm = col*rm; 200807171c5SHong Zhang for (i=0; i<rm; i++) { /* over rows of R in those columns */ 201807171c5SHong Zhang r1 = r2 = r3 = r4 = 0.0; 202807171c5SHong Zhang n = r->i[i+1] - r->i[i]; 203807171c5SHong Zhang rj = r->j + r->i[i]; 204807171c5SHong Zhang ra = r->a + r->i[i]; 205807171c5SHong Zhang for (j=0; j<n; j++) { 206807171c5SHong Zhang r1 += (*ra)*c[*rj]; 207807171c5SHong Zhang r2 += (*ra)*c2[*rj]; 208807171c5SHong Zhang r3 += (*ra)*c3[*rj]; 209807171c5SHong Zhang r4 += (*ra++)*c4[*rj++]; 210807171c5SHong Zhang } 211807171c5SHong Zhang d[colrm + i] = r1; 212807171c5SHong Zhang d[colrm + rm + i] = r2; 213807171c5SHong Zhang d[colrm + rm2 + i] = r3; 214807171c5SHong Zhang d[colrm + rm3 + i] = r4; 215807171c5SHong Zhang } 216807171c5SHong Zhang } 217807171c5SHong Zhang for (; col<cn; col++) { /* over extra columns of C */ 218807171c5SHong Zhang for (i=0; i<am; i++) { /* over rows of A in those columns */ 219807171c5SHong Zhang r1 = 0.0; 220807171c5SHong Zhang n = a->i[i+1] - a->i[i]; 221807171c5SHong Zhang aj = a->j + a->i[i]; 222807171c5SHong Zhang aa = a->a + a->i[i]; 223807171c5SHong Zhang for (j=0; j<n; j++) { 224807171c5SHong Zhang r1 += (*aa++)*b1[*aj++]; 225807171c5SHong Zhang } 226807171c5SHong Zhang c[i] = r1; 227807171c5SHong Zhang } 228807171c5SHong Zhang b1 += bm; 229807171c5SHong Zhang 230807171c5SHong Zhang for (i=0; i<rm; i++) { /* over rows of R in those columns */ 231807171c5SHong Zhang r1 = 0.0; 232807171c5SHong Zhang n = r->i[i+1] - r->i[i]; 233807171c5SHong Zhang rj = r->j + r->i[i]; 234807171c5SHong Zhang ra = r->a + r->i[i]; 235807171c5SHong Zhang for (j=0; j<n; j++) { 236807171c5SHong Zhang r1 += (*ra++)*c[*rj++]; 237807171c5SHong Zhang } 238807171c5SHong Zhang d[col*rm + i] = r1; 239807171c5SHong Zhang } 240807171c5SHong Zhang } 241807171c5SHong Zhang ierr = PetscLogFlops(cn*2.0*(a->nz + r->nz));CHKERRQ(ierr); 242807171c5SHong Zhang 2438c778c55SBarry Smith ierr = MatDenseRestoreArray(B,&b);CHKERRQ(ierr); 2448c778c55SBarry Smith ierr = MatDenseRestoreArray(RAB,&d);CHKERRQ(ierr); 245807171c5SHong Zhang ierr = MatAssemblyBegin(RAB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 246807171c5SHong Zhang ierr = MatAssemblyEnd(RAB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 247807171c5SHong Zhang PetscFunctionReturn(0); 248807171c5SHong Zhang } 249807171c5SHong Zhang 250807171c5SHong Zhang #undef __FUNCT__ 251807171c5SHong Zhang #define __FUNCT__ "MatRARtNumeric_SeqAIJ_SeqAIJ" 252807171c5SHong Zhang PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ(Mat A,Mat R,Mat C) 253807171c5SHong Zhang { 254807171c5SHong Zhang PetscErrorCode ierr; 255257c235dSHong Zhang Mat_SeqAIJ *c = (Mat_SeqAIJ*)C->data; 256257c235dSHong Zhang Mat_RARt *rart=c->rart; 257807171c5SHong Zhang MatTransposeColoring matcoloring; 2588d0a38e4SHong Zhang Mat Rt,RARt; 259807171c5SHong Zhang PetscLogDouble Mult_sp_den=0.0,app1=0.0,app2=0.0,t0,tf; 260807171c5SHong Zhang 261807171c5SHong Zhang PetscFunctionBegin; 262807171c5SHong Zhang /* Get dense Rt by Apply MatTransposeColoring to R */ 263807171c5SHong Zhang matcoloring = rart->matcoloring; 264807171c5SHong Zhang Rt = rart->Rt; 2652205254eSKarl Rupp 2668563dfccSBarry Smith ierr = PetscTime(&t0);CHKERRQ(ierr); 267807171c5SHong Zhang ierr = MatTransColoringApplySpToDen(matcoloring,R,Rt);CHKERRQ(ierr); 2688563dfccSBarry Smith ierr = PetscTime(&tf);CHKERRQ(ierr); 269807171c5SHong Zhang app1 += tf - t0; 270807171c5SHong Zhang 27150647e95SHong Zhang /* Get dense RARt = R*A*Rt -- dominates! */ 2728563dfccSBarry Smith ierr = PetscTime(&t0);CHKERRQ(ierr); 273807171c5SHong Zhang RARt = rart->RARt; 274807171c5SHong Zhang ierr = MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense(R,A,Rt,RARt,rart->work);CHKERRQ(ierr); 2758563dfccSBarry Smith ierr = PetscTime(&tf);CHKERRQ(ierr); 2762205254eSKarl Rupp 277807171c5SHong Zhang Mult_sp_den += tf - t0; 278807171c5SHong Zhang 279807171c5SHong Zhang /* Recover C from C_dense */ 2808563dfccSBarry Smith ierr = PetscTime(&t0);CHKERRQ(ierr); 281807171c5SHong Zhang ierr = MatTransColoringApplyDenToSp(matcoloring,RARt,C);CHKERRQ(ierr); 2828563dfccSBarry Smith ierr = PetscTime(&tf);CHKERRQ(ierr); 2832205254eSKarl Rupp 284807171c5SHong Zhang app2 += tf - t0; 285807171c5SHong Zhang 286807171c5SHong Zhang #if defined(PETSC_USE_INFO) 287a2ea699eSBarry Smith ierr = PetscInfo4(C,"Num = ColorApp %g + %g + Mult_sp_den %g = %g\n",app1,app2,Mult_sp_den,app1+app2+Mult_sp_den);CHKERRQ(ierr); 288807171c5SHong Zhang #endif 289807171c5SHong Zhang PetscFunctionReturn(0); 290807171c5SHong Zhang } 291807171c5SHong Zhang 292807171c5SHong Zhang #undef __FUNCT__ 293807171c5SHong Zhang #define __FUNCT__ "MatRARt_SeqAIJ_SeqAIJ" 294807171c5SHong Zhang PetscErrorCode MatRARt_SeqAIJ_SeqAIJ(Mat A,Mat R,MatReuse scall,PetscReal fill,Mat *C) 295807171c5SHong Zhang { 296807171c5SHong Zhang PetscErrorCode ierr; 297*4fa85fdeSHong Zhang PetscBool usecoloring = PETSC_FALSE; 298807171c5SHong Zhang 299807171c5SHong Zhang PetscFunctionBegin; 30022e94b5dSHong Zhang ierr = PetscOptionsGetBool(NULL,"-matrart_color",&usecoloring,NULL);CHKERRQ(ierr); 30122e94b5dSHong Zhang 302*4fa85fdeSHong Zhang if (!usecoloring) { /* Rt=R^T, C=R*A*Rt - avoid inefficient sparse inner products */ 30322e94b5dSHong Zhang Mat Rt; 304*4fa85fdeSHong Zhang Mat_SeqAIJ *c; 305*4fa85fdeSHong Zhang Mat_RARt *rart; 306*4fa85fdeSHong Zhang 307*4fa85fdeSHong Zhang if (scall == MAT_INITIAL_MATRIX) { 308*4fa85fdeSHong Zhang ierr = MatTranspose(R,MAT_INITIAL_MATRIX,&Rt);CHKERRQ(ierr); 30922e94b5dSHong Zhang ierr = MatMatMatMult(R,A,Rt,scall,fill,C);CHKERRQ(ierr); 310*4fa85fdeSHong Zhang 311*4fa85fdeSHong Zhang ierr = PetscNew(Mat_RARt,&rart);CHKERRQ(ierr); 312*4fa85fdeSHong Zhang rart->Rt = Rt; 313*4fa85fdeSHong Zhang c = (Mat_SeqAIJ*)(*C)->data; 314*4fa85fdeSHong Zhang c->rart = rart; 315*4fa85fdeSHong Zhang rart->destroy = (*C)->ops->destroy; 316*4fa85fdeSHong Zhang (*C)->ops->destroy = MatDestroy_SeqAIJ_RARt; 317*4fa85fdeSHong Zhang } else { 318*4fa85fdeSHong Zhang c = (Mat_SeqAIJ*)(*C)->data; 319*4fa85fdeSHong Zhang rart = c->rart; 320*4fa85fdeSHong Zhang Rt = rart->Rt; 321*4fa85fdeSHong Zhang ierr = MatTranspose(R,MAT_REUSE_MATRIX,&Rt);CHKERRQ(ierr); 322*4fa85fdeSHong Zhang ierr = MatMatMatMult(R,A,Rt,scall,fill,C);CHKERRQ(ierr); 323*4fa85fdeSHong Zhang } 32422e94b5dSHong Zhang PetscFunctionReturn(0); 32522e94b5dSHong Zhang } 32622e94b5dSHong Zhang 32722e94b5dSHong Zhang /* use coloring */ 32836104f73SHong Zhang PetscBool newimpl=PETSC_FALSE; 32936104f73SHong Zhang ierr = PetscOptionsGetBool(NULL,"-matrart_new",&newimpl,NULL);CHKERRQ(ierr); 33036104f73SHong Zhang if (newimpl) { 33136104f73SHong Zhang Mat ARt,RARt; 33236104f73SHong Zhang printf(" MatRARt_new ....\n"); 33336104f73SHong Zhang if (scall == MAT_INITIAL_MATRIX) { 33436104f73SHong Zhang ierr = PetscLogEventBegin(MAT_RARtSymbolic,A,R,0,0);CHKERRQ(ierr); 33536104f73SHong Zhang ierr = MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(A,R,fill,&ARt);CHKERRQ(ierr); 33636104f73SHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(R,ARt,fill,&RARt);CHKERRQ(ierr); 33736104f73SHong Zhang ierr = PetscLogEventEnd(MAT_RARtSymbolic,A,R,0,0);CHKERRQ(ierr); 33836104f73SHong Zhang } 33936104f73SHong Zhang ierr = PetscLogEventBegin(MAT_RARtNumeric,A,R,0,0);CHKERRQ(ierr); 34036104f73SHong Zhang ierr = MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(A,R,ARt);CHKERRQ(ierr); 34136104f73SHong Zhang ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ(R,ARt,RARt);CHKERRQ(ierr); 34236104f73SHong Zhang *C = RARt; 34336104f73SHong Zhang ierr = PetscLogEventEnd(MAT_RARtNumeric,A,R,0,0);CHKERRQ(ierr); 34436104f73SHong Zhang 34536104f73SHong Zhang ierr = MatDestroy(&ARt);CHKERRQ(ierr); 34636104f73SHong Zhang PetscFunctionReturn(0); 34736104f73SHong Zhang } 348807171c5SHong Zhang if (scall == MAT_INITIAL_MATRIX) { 34950647e95SHong Zhang ierr = PetscLogEventBegin(MAT_RARtSymbolic,A,R,0,0);CHKERRQ(ierr); 350807171c5SHong Zhang ierr = MatRARtSymbolic_SeqAIJ_SeqAIJ(A,R,fill,C);CHKERRQ(ierr); 35150647e95SHong Zhang ierr = PetscLogEventEnd(MAT_RARtSymbolic,A,R,0,0);CHKERRQ(ierr); 352807171c5SHong Zhang } 35350647e95SHong Zhang ierr = PetscLogEventBegin(MAT_RARtNumeric,A,R,0,0);CHKERRQ(ierr); 354807171c5SHong Zhang ierr = MatRARtNumeric_SeqAIJ_SeqAIJ(A,R,*C);CHKERRQ(ierr); 35550647e95SHong Zhang ierr = PetscLogEventEnd(MAT_RARtNumeric,A,R,0,0);CHKERRQ(ierr); 356807171c5SHong Zhang PetscFunctionReturn(0); 357807171c5SHong Zhang } 358