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__ 13807171c5SHong Zhang #define __FUNCT__ "PetscContainerDestroy_Mat_RARt" 14807171c5SHong Zhang PetscErrorCode PetscContainerDestroy_Mat_RARt(void *ptr) 15807171c5SHong Zhang { 16807171c5SHong Zhang PetscErrorCode ierr; 17807171c5SHong Zhang Mat_RARt *rart=(Mat_RARt*)ptr; 18807171c5SHong Zhang 19807171c5SHong Zhang PetscFunctionBegin; 20807171c5SHong Zhang ierr = MatTransposeColoringDestroy(&rart->matcoloring);CHKERRQ(ierr); 21807171c5SHong Zhang ierr = MatDestroy(&rart->Rt);CHKERRQ(ierr); 22807171c5SHong Zhang ierr = MatDestroy(&rart->RARt);CHKERRQ(ierr); 23807171c5SHong Zhang ierr = PetscFree(rart->work);CHKERRQ(ierr); 24807171c5SHong Zhang ierr = PetscFree(rart);CHKERRQ(ierr); 25807171c5SHong Zhang PetscFunctionReturn(0); 26807171c5SHong Zhang } 27807171c5SHong Zhang 28807171c5SHong Zhang #undef __FUNCT__ 29807171c5SHong Zhang #define __FUNCT__ "MatDestroy_SeqAIJ_RARt" 30807171c5SHong Zhang PetscErrorCode MatDestroy_SeqAIJ_RARt(Mat A) 31807171c5SHong Zhang { 32807171c5SHong Zhang PetscErrorCode ierr; 33807171c5SHong Zhang PetscContainer container; 340298fd71SBarry Smith Mat_RARt *rart=NULL; 35807171c5SHong Zhang 36807171c5SHong Zhang PetscFunctionBegin; 37807171c5SHong Zhang ierr = PetscObjectQuery((PetscObject)A,"Mat_RARt",(PetscObject*)&container);CHKERRQ(ierr); 38807171c5SHong Zhang if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit"); 39807171c5SHong Zhang ierr = PetscContainerGetPointer(container,(void**)&rart);CHKERRQ(ierr); 40807171c5SHong Zhang A->ops->destroy = rart->destroy; 41807171c5SHong Zhang if (A->ops->destroy) { 42807171c5SHong Zhang ierr = (*A->ops->destroy)(A);CHKERRQ(ierr); 43807171c5SHong Zhang } 44807171c5SHong Zhang ierr = PetscObjectCompose((PetscObject)A,"Mat_RARt",0);CHKERRQ(ierr); 45807171c5SHong Zhang PetscFunctionReturn(0); 46807171c5SHong Zhang } 47807171c5SHong Zhang 48807171c5SHong Zhang #undef __FUNCT__ 49807171c5SHong Zhang #define __FUNCT__ "MatRARtSymbolic_SeqAIJ_SeqAIJ" 50807171c5SHong Zhang PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat R,PetscReal fill,Mat *C) 51807171c5SHong Zhang { 52807171c5SHong Zhang PetscErrorCode ierr; 53807171c5SHong Zhang Mat P; 54807171c5SHong Zhang PetscInt *rti,*rtj; 55807171c5SHong Zhang Mat_RARt *rart; 56807171c5SHong Zhang PetscContainer container; 57807171c5SHong Zhang MatTransposeColoring matcoloring; 58807171c5SHong Zhang ISColoring iscoloring; 598d0a38e4SHong Zhang Mat Rt_dense,RARt_dense; 60807171c5SHong Zhang PetscLogDouble GColor=0.0,MCCreate=0.0,MDenCreate=0.0,t0,tf,etime=0.0; 61807171c5SHong Zhang Mat_SeqAIJ *c; 62807171c5SHong Zhang 63807171c5SHong Zhang PetscFunctionBegin; 648563dfccSBarry Smith ierr = PetscTime(&t0);CHKERRQ(ierr); 65807171c5SHong Zhang /* create symbolic P=Rt */ 66807171c5SHong Zhang ierr = MatGetSymbolicTranspose_SeqAIJ(R,&rti,&rtj);CHKERRQ(ierr); 670298fd71SBarry Smith ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,R->cmap->n,R->rmap->n,rti,rtj,NULL,&P);CHKERRQ(ierr); 68807171c5SHong Zhang 69807171c5SHong Zhang /* get symbolic C=Pt*A*P */ 70807171c5SHong Zhang ierr = MatPtAPSymbolic_SeqAIJ_SeqAIJ(A,P,fill,C);CHKERRQ(ierr); 71fd8e3dafSHong Zhang (*C)->rmap->bs = R->rmap->bs; 72fd8e3dafSHong Zhang (*C)->cmap->bs = R->rmap->bs; 73807171c5SHong Zhang 74807171c5SHong Zhang /* create a supporting struct */ 75807171c5SHong Zhang ierr = PetscNew(Mat_RARt,&rart);CHKERRQ(ierr); 76807171c5SHong Zhang 77807171c5SHong Zhang /* attach the supporting struct to C */ 78807171c5SHong Zhang ierr = PetscContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr); 79807171c5SHong Zhang ierr = PetscContainerSetPointer(container,rart);CHKERRQ(ierr); 80807171c5SHong Zhang ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_Mat_RARt);CHKERRQ(ierr); 81807171c5SHong Zhang ierr = PetscObjectCompose((PetscObject)(*C),"Mat_RARt",(PetscObject)container);CHKERRQ(ierr); 82807171c5SHong Zhang ierr = PetscContainerDestroy(&container);CHKERRQ(ierr); 838563dfccSBarry Smith ierr = PetscTime(&tf);CHKERRQ(ierr); 84807171c5SHong Zhang etime += tf - t0; 85807171c5SHong Zhang 8622e94b5dSHong Zhang /* ------ Use coloring ---------- */ 87807171c5SHong Zhang /* Create MatTransposeColoring from symbolic C=R*A*R^T */ 88807171c5SHong Zhang c = (Mat_SeqAIJ*)(*C)->data; 898563dfccSBarry Smith ierr = PetscTime(&t0);CHKERRQ(ierr); 90807171c5SHong Zhang ierr = MatGetColoring(*C,MATCOLORINGLF,&iscoloring);CHKERRQ(ierr); 918563dfccSBarry Smith ierr = PetscTime(&tf);CHKERRQ(ierr); 92807171c5SHong Zhang GColor += tf - t0; 93807171c5SHong Zhang 948563dfccSBarry Smith ierr = PetscTime(&t0);CHKERRQ(ierr); 95807171c5SHong Zhang ierr = MatTransposeColoringCreate(*C,iscoloring,&matcoloring);CHKERRQ(ierr); 962205254eSKarl Rupp 97807171c5SHong Zhang rart->matcoloring = matcoloring; 982205254eSKarl Rupp 99807171c5SHong Zhang ierr = ISColoringDestroy(&iscoloring);CHKERRQ(ierr); 1008563dfccSBarry Smith ierr = PetscTime(&tf);CHKERRQ(ierr); 101807171c5SHong Zhang MCCreate += tf - t0; 102807171c5SHong Zhang 1038563dfccSBarry Smith ierr = PetscTime(&t0);CHKERRQ(ierr); 104807171c5SHong Zhang /* Create Rt_dense */ 105807171c5SHong Zhang ierr = MatCreate(PETSC_COMM_SELF,&Rt_dense);CHKERRQ(ierr); 106807171c5SHong Zhang ierr = MatSetSizes(Rt_dense,A->cmap->n,matcoloring->ncolors,A->cmap->n,matcoloring->ncolors);CHKERRQ(ierr); 107807171c5SHong Zhang ierr = MatSetType(Rt_dense,MATSEQDENSE);CHKERRQ(ierr); 1080298fd71SBarry Smith ierr = MatSeqDenseSetPreallocation(Rt_dense,NULL);CHKERRQ(ierr); 1092205254eSKarl Rupp 110807171c5SHong Zhang Rt_dense->assembled = PETSC_TRUE; 111807171c5SHong Zhang rart->Rt = Rt_dense; 112807171c5SHong Zhang 113807171c5SHong Zhang /* Create RARt_dense = R*A*Rt_dense */ 114807171c5SHong Zhang ierr = MatCreate(PETSC_COMM_SELF,&RARt_dense);CHKERRQ(ierr); 115807171c5SHong Zhang ierr = MatSetSizes(RARt_dense,(*C)->rmap->n,matcoloring->ncolors,(*C)->rmap->n,matcoloring->ncolors);CHKERRQ(ierr); 116807171c5SHong Zhang ierr = MatSetType(RARt_dense,MATSEQDENSE);CHKERRQ(ierr); 1170298fd71SBarry Smith ierr = MatSeqDenseSetPreallocation(RARt_dense,NULL);CHKERRQ(ierr); 1182205254eSKarl Rupp 119807171c5SHong Zhang rart->RARt = RARt_dense; 120807171c5SHong Zhang 121807171c5SHong Zhang /* Allocate work array to store columns of A*R^T used in MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense() */ 122807171c5SHong Zhang ierr = PetscMalloc(A->rmap->n*4*sizeof(PetscScalar),&rart->work);CHKERRQ(ierr); 123807171c5SHong Zhang 1248563dfccSBarry Smith ierr = PetscTime(&tf);CHKERRQ(ierr); 125807171c5SHong Zhang MDenCreate += tf - t0; 126807171c5SHong Zhang 127807171c5SHong Zhang rart->destroy = (*C)->ops->destroy; 128807171c5SHong Zhang (*C)->ops->destroy = MatDestroy_SeqAIJ_RARt; 129807171c5SHong Zhang 130807171c5SHong Zhang /* clean up */ 131807171c5SHong Zhang ierr = MatRestoreSymbolicTranspose_SeqAIJ(R,&rti,&rtj);CHKERRQ(ierr); 132807171c5SHong Zhang ierr = MatDestroy(&P);CHKERRQ(ierr); 133807171c5SHong Zhang 134807171c5SHong Zhang #if defined(PETSC_USE_INFO) 1351ce71dffSSatish Balay { 136807171c5SHong Zhang PetscReal density= (PetscReal)(c->nz)/(RARt_dense->rmap->n*RARt_dense->cmap->n); 137a2ea699eSBarry 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); 138a2ea699eSBarry 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); 1391ce71dffSSatish Balay } 140807171c5SHong Zhang #endif 141807171c5SHong Zhang PetscFunctionReturn(0); 142807171c5SHong Zhang } 143807171c5SHong Zhang 144807171c5SHong Zhang /* 145807171c5SHong Zhang RAB = R * A * B, R and A in seqaij format, B in dense format; 146807171c5SHong Zhang */ 147807171c5SHong Zhang #undef __FUNCT__ 148807171c5SHong Zhang #define __FUNCT__ "MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense" 149807171c5SHong Zhang PetscErrorCode MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense(Mat R,Mat A,Mat B,Mat RAB,PetscScalar *work) 150807171c5SHong Zhang { 151807171c5SHong Zhang Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*r=(Mat_SeqAIJ*)R->data; 152807171c5SHong Zhang PetscErrorCode ierr; 153807171c5SHong Zhang PetscScalar *b,r1,r2,r3,r4,*b1,*b2,*b3,*b4; 154807171c5SHong Zhang MatScalar *aa,*ra; 155807171c5SHong Zhang PetscInt cn =B->cmap->n,bm=B->rmap->n,col,i,j,n,*ai=a->i,*aj,am=A->rmap->n; 156807171c5SHong Zhang PetscInt am2=2*am,am3=3*am,bm4=4*bm; 157807171c5SHong Zhang PetscScalar *d,*c,*c2,*c3,*c4; 158807171c5SHong Zhang PetscInt *rj,rm=R->rmap->n,dm=RAB->rmap->n,dn=RAB->cmap->n; 159807171c5SHong Zhang PetscInt rm2=2*rm,rm3=3*rm,colrm; 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 1688c778c55SBarry Smith ierr = MatDenseGetArray(B,&b);CHKERRQ(ierr); 1698c778c55SBarry Smith ierr = MatDenseGetArray(RAB,&d);CHKERRQ(ierr); 170807171c5SHong Zhang b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm; 171807171c5SHong Zhang c = work; c2 = c + am; c3 = c2 + am; c4 = c3 + am; 172807171c5SHong Zhang for (col=0; col<cn-4; col += 4) { /* over columns of C */ 173807171c5SHong Zhang for (i=0; i<am; i++) { /* over rows of A in those columns */ 174807171c5SHong Zhang r1 = r2 = r3 = r4 = 0.0; 175807171c5SHong Zhang n = ai[i+1] - ai[i]; 176807171c5SHong Zhang aj = a->j + ai[i]; 177807171c5SHong Zhang aa = a->a + ai[i]; 178807171c5SHong Zhang for (j=0; j<n; j++) { 179807171c5SHong Zhang r1 += (*aa)*b1[*aj]; 180807171c5SHong Zhang r2 += (*aa)*b2[*aj]; 181807171c5SHong Zhang r3 += (*aa)*b3[*aj]; 182807171c5SHong Zhang r4 += (*aa++)*b4[*aj++]; 183807171c5SHong Zhang } 184807171c5SHong Zhang c[i] = r1; 185807171c5SHong Zhang c[am + i] = r2; 186807171c5SHong Zhang c[am2 + i] = r3; 187807171c5SHong Zhang c[am3 + i] = r4; 188807171c5SHong Zhang } 189807171c5SHong Zhang b1 += bm4; 190807171c5SHong Zhang b2 += bm4; 191807171c5SHong Zhang b3 += bm4; 192807171c5SHong Zhang b4 += bm4; 193807171c5SHong Zhang 194807171c5SHong Zhang /* RAB[:,col] = R*C[:,col] */ 195807171c5SHong Zhang colrm = col*rm; 196807171c5SHong Zhang for (i=0; i<rm; i++) { /* over rows of R in those columns */ 197807171c5SHong Zhang r1 = r2 = r3 = r4 = 0.0; 198807171c5SHong Zhang n = r->i[i+1] - r->i[i]; 199807171c5SHong Zhang rj = r->j + r->i[i]; 200807171c5SHong Zhang ra = r->a + r->i[i]; 201807171c5SHong Zhang for (j=0; j<n; j++) { 202807171c5SHong Zhang r1 += (*ra)*c[*rj]; 203807171c5SHong Zhang r2 += (*ra)*c2[*rj]; 204807171c5SHong Zhang r3 += (*ra)*c3[*rj]; 205807171c5SHong Zhang r4 += (*ra++)*c4[*rj++]; 206807171c5SHong Zhang } 207807171c5SHong Zhang d[colrm + i] = r1; 208807171c5SHong Zhang d[colrm + rm + i] = r2; 209807171c5SHong Zhang d[colrm + rm2 + i] = r3; 210807171c5SHong Zhang d[colrm + rm3 + i] = r4; 211807171c5SHong Zhang } 212807171c5SHong Zhang } 213807171c5SHong Zhang for (; col<cn; col++) { /* over extra columns of C */ 214807171c5SHong Zhang for (i=0; i<am; i++) { /* over rows of A in those columns */ 215807171c5SHong Zhang r1 = 0.0; 216807171c5SHong Zhang n = a->i[i+1] - a->i[i]; 217807171c5SHong Zhang aj = a->j + a->i[i]; 218807171c5SHong Zhang aa = a->a + a->i[i]; 219807171c5SHong Zhang for (j=0; j<n; j++) { 220807171c5SHong Zhang r1 += (*aa++)*b1[*aj++]; 221807171c5SHong Zhang } 222807171c5SHong Zhang c[i] = r1; 223807171c5SHong Zhang } 224807171c5SHong Zhang b1 += bm; 225807171c5SHong Zhang 226807171c5SHong Zhang for (i=0; i<rm; i++) { /* over rows of R in those columns */ 227807171c5SHong Zhang r1 = 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 } 234807171c5SHong Zhang d[col*rm + i] = r1; 235807171c5SHong Zhang } 236807171c5SHong Zhang } 237807171c5SHong Zhang ierr = PetscLogFlops(cn*2.0*(a->nz + r->nz));CHKERRQ(ierr); 238807171c5SHong Zhang 2398c778c55SBarry Smith ierr = MatDenseRestoreArray(B,&b);CHKERRQ(ierr); 2408c778c55SBarry Smith ierr = MatDenseRestoreArray(RAB,&d);CHKERRQ(ierr); 241807171c5SHong Zhang ierr = MatAssemblyBegin(RAB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 242807171c5SHong Zhang ierr = MatAssemblyEnd(RAB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 243807171c5SHong Zhang PetscFunctionReturn(0); 244807171c5SHong Zhang } 245807171c5SHong Zhang 246807171c5SHong Zhang #undef __FUNCT__ 247807171c5SHong Zhang #define __FUNCT__ "MatRARtNumeric_SeqAIJ_SeqAIJ" 248807171c5SHong Zhang PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ(Mat A,Mat R,Mat C) 249807171c5SHong Zhang { 250807171c5SHong Zhang PetscErrorCode ierr; 251807171c5SHong Zhang Mat_RARt *rart; 252807171c5SHong Zhang PetscContainer container; 253807171c5SHong Zhang MatTransposeColoring matcoloring; 2548d0a38e4SHong Zhang Mat Rt,RARt; 255807171c5SHong Zhang PetscLogDouble Mult_sp_den=0.0,app1=0.0,app2=0.0,t0,tf; 256807171c5SHong Zhang 257807171c5SHong Zhang PetscFunctionBegin; 258807171c5SHong Zhang ierr = PetscObjectQuery((PetscObject)C,"Mat_RARt",(PetscObject*)&container);CHKERRQ(ierr); 259807171c5SHong Zhang if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit"); 260807171c5SHong Zhang ierr = PetscContainerGetPointer(container,(void**)&rart);CHKERRQ(ierr); 261807171c5SHong Zhang 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; 29722e94b5dSHong Zhang PetscBool usecoloring = PETSC_TRUE; 298807171c5SHong Zhang 299807171c5SHong Zhang PetscFunctionBegin; 30022e94b5dSHong Zhang ierr = PetscOptionsGetBool(NULL,"-matrart_color",&usecoloring,NULL);CHKERRQ(ierr); 30122e94b5dSHong Zhang 30222e94b5dSHong Zhang if (!usecoloring) { 30322e94b5dSHong Zhang Mat Rt; 30422e94b5dSHong Zhang ierr = MatTranspose(R,MAT_INITIAL_MATRIX,&Rt);CHKERRQ(ierr); /* replace MAT_INITIAL_MATRIX with scall if !usecoloring is better */ 30522e94b5dSHong Zhang ierr = MatMatMatMult(R,A,Rt,scall,fill,C);CHKERRQ(ierr); 30622e94b5dSHong Zhang ierr = MatDestroy(&Rt);CHKERRQ(ierr); 30722e94b5dSHong Zhang PetscFunctionReturn(0); 30822e94b5dSHong Zhang } 30922e94b5dSHong Zhang 31022e94b5dSHong Zhang /* use coloring */ 311*36104f73SHong Zhang PetscBool newimpl=PETSC_FALSE; 312*36104f73SHong Zhang ierr = PetscOptionsGetBool(NULL,"-matrart_new",&newimpl,NULL);CHKERRQ(ierr); 313*36104f73SHong Zhang if (newimpl) { 314*36104f73SHong Zhang Mat ARt,RARt; 315*36104f73SHong Zhang printf(" MatRARt_new ....\n"); 316*36104f73SHong Zhang if (scall == MAT_INITIAL_MATRIX) { 317*36104f73SHong Zhang ierr = PetscLogEventBegin(MAT_RARtSymbolic,A,R,0,0);CHKERRQ(ierr); 318*36104f73SHong Zhang ierr = MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(A,R,fill,&ARt);CHKERRQ(ierr); 319*36104f73SHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(R,ARt,fill,&RARt);CHKERRQ(ierr); 320*36104f73SHong Zhang ierr = PetscLogEventEnd(MAT_RARtSymbolic,A,R,0,0);CHKERRQ(ierr); 321*36104f73SHong Zhang } 322*36104f73SHong Zhang ierr = PetscLogEventBegin(MAT_RARtNumeric,A,R,0,0);CHKERRQ(ierr); 323*36104f73SHong Zhang ierr = MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(A,R,ARt);CHKERRQ(ierr); 324*36104f73SHong Zhang ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ(R,ARt,RARt);CHKERRQ(ierr); 325*36104f73SHong Zhang *C = RARt; 326*36104f73SHong Zhang ierr = PetscLogEventEnd(MAT_RARtNumeric,A,R,0,0);CHKERRQ(ierr); 327*36104f73SHong Zhang 328*36104f73SHong Zhang ierr = MatDestroy(&ARt);CHKERRQ(ierr); 329*36104f73SHong Zhang PetscFunctionReturn(0); 330*36104f73SHong Zhang } 331807171c5SHong Zhang if (scall == MAT_INITIAL_MATRIX) { 33250647e95SHong Zhang ierr = PetscLogEventBegin(MAT_RARtSymbolic,A,R,0,0);CHKERRQ(ierr); 333807171c5SHong Zhang ierr = MatRARtSymbolic_SeqAIJ_SeqAIJ(A,R,fill,C);CHKERRQ(ierr); 33450647e95SHong Zhang ierr = PetscLogEventEnd(MAT_RARtSymbolic,A,R,0,0);CHKERRQ(ierr); 335807171c5SHong Zhang } 33650647e95SHong Zhang ierr = PetscLogEventBegin(MAT_RARtNumeric,A,R,0,0);CHKERRQ(ierr); 337807171c5SHong Zhang ierr = MatRARtNumeric_SeqAIJ_SeqAIJ(A,R,*C);CHKERRQ(ierr); 33850647e95SHong Zhang ierr = PetscLogEventEnd(MAT_RARtNumeric,A,R,0,0);CHKERRQ(ierr); 339807171c5SHong Zhang PetscFunctionReturn(0); 340807171c5SHong Zhang } 341