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 11257c235dSHong Zhang PetscErrorCode MatDestroy_SeqAIJ_RARt(Mat A) 12807171c5SHong Zhang { 13807171c5SHong Zhang PetscErrorCode ierr; 14257c235dSHong Zhang Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 15257c235dSHong Zhang Mat_RARt *rart = a->rart; 16807171c5SHong Zhang 17807171c5SHong Zhang PetscFunctionBegin; 18807171c5SHong Zhang ierr = MatTransposeColoringDestroy(&rart->matcoloring);CHKERRQ(ierr); 19807171c5SHong Zhang ierr = MatDestroy(&rart->Rt);CHKERRQ(ierr); 20807171c5SHong Zhang ierr = MatDestroy(&rart->RARt);CHKERRQ(ierr); 213b1b9624SHong Zhang ierr = MatDestroy(&rart->ARt);CHKERRQ(ierr); 22807171c5SHong Zhang ierr = PetscFree(rart->work);CHKERRQ(ierr); 23807171c5SHong Zhang 24807171c5SHong Zhang A->ops->destroy = rart->destroy; 25807171c5SHong Zhang if (A->ops->destroy) { 26807171c5SHong Zhang ierr = (*A->ops->destroy)(A);CHKERRQ(ierr); 27807171c5SHong Zhang } 28257c235dSHong Zhang ierr = PetscFree(rart);CHKERRQ(ierr); 29807171c5SHong Zhang PetscFunctionReturn(0); 30807171c5SHong Zhang } 31807171c5SHong Zhang 3255bea0ebSHong Zhang PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ_colorrart(Mat A,Mat R,PetscReal fill,Mat *C) 33807171c5SHong Zhang { 34807171c5SHong Zhang PetscErrorCode ierr; 35807171c5SHong Zhang Mat P; 36807171c5SHong Zhang PetscInt *rti,*rtj; 37807171c5SHong Zhang Mat_RARt *rart; 38335efc43SPeter Brune MatColoring coloring; 39807171c5SHong Zhang MatTransposeColoring matcoloring; 40807171c5SHong Zhang ISColoring iscoloring; 418d0a38e4SHong Zhang Mat Rt_dense,RARt_dense; 42807171c5SHong Zhang Mat_SeqAIJ *c; 43807171c5SHong Zhang 44807171c5SHong Zhang PetscFunctionBegin; 45807171c5SHong Zhang /* create symbolic P=Rt */ 46807171c5SHong Zhang ierr = MatGetSymbolicTranspose_SeqAIJ(R,&rti,&rtj);CHKERRQ(ierr); 470298fd71SBarry Smith ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,R->cmap->n,R->rmap->n,rti,rtj,NULL,&P);CHKERRQ(ierr); 48807171c5SHong Zhang 49807171c5SHong Zhang /* get symbolic C=Pt*A*P */ 504a1b09b7SHong Zhang ierr = MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy(A,P,fill,C);CHKERRQ(ierr); 5133d57670SJed Brown ierr = MatSetBlockSizes(*C,PetscAbs(R->rmap->bs),PetscAbs(R->rmap->bs));CHKERRQ(ierr); 5255bea0ebSHong Zhang (*C)->ops->rartnumeric = MatRARtNumeric_SeqAIJ_SeqAIJ_colorrart; 53807171c5SHong Zhang 54807171c5SHong Zhang /* create a supporting struct */ 55b00a9115SJed Brown ierr = PetscNew(&rart);CHKERRQ(ierr); 56257c235dSHong Zhang c = (Mat_SeqAIJ*)(*C)->data; 57257c235dSHong Zhang c->rart = rart; 58807171c5SHong Zhang 5922e94b5dSHong Zhang /* ------ Use coloring ---------- */ 604d478ae7SHong Zhang /* inode causes memory problem, don't know why */ 614d478ae7SHong Zhang if (c->inode.use) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MAT_USE_INODES is not supported. Use '-mat_no_inode'"); 624d478ae7SHong Zhang 63807171c5SHong Zhang /* Create MatTransposeColoring from symbolic C=R*A*R^T */ 64335efc43SPeter Brune ierr = MatColoringCreate(*C,&coloring);CHKERRQ(ierr); 65335efc43SPeter Brune ierr = MatColoringSetDistance(coloring,2);CHKERRQ(ierr); 66335efc43SPeter Brune ierr = MatColoringSetType(coloring,MATCOLORINGSL);CHKERRQ(ierr); 67335efc43SPeter Brune ierr = MatColoringSetFromOptions(coloring);CHKERRQ(ierr); 68335efc43SPeter Brune ierr = MatColoringApply(coloring,&iscoloring);CHKERRQ(ierr); 69335efc43SPeter Brune ierr = MatColoringDestroy(&coloring);CHKERRQ(ierr); 70807171c5SHong Zhang ierr = MatTransposeColoringCreate(*C,iscoloring,&matcoloring);CHKERRQ(ierr); 712205254eSKarl Rupp 72807171c5SHong Zhang rart->matcoloring = matcoloring; 73807171c5SHong Zhang ierr = ISColoringDestroy(&iscoloring);CHKERRQ(ierr); 743b1b9624SHong Zhang 75807171c5SHong Zhang /* Create Rt_dense */ 76807171c5SHong Zhang ierr = MatCreate(PETSC_COMM_SELF,&Rt_dense);CHKERRQ(ierr); 77807171c5SHong Zhang ierr = MatSetSizes(Rt_dense,A->cmap->n,matcoloring->ncolors,A->cmap->n,matcoloring->ncolors);CHKERRQ(ierr); 78807171c5SHong Zhang ierr = MatSetType(Rt_dense,MATSEQDENSE);CHKERRQ(ierr); 790298fd71SBarry Smith ierr = MatSeqDenseSetPreallocation(Rt_dense,NULL);CHKERRQ(ierr); 802205254eSKarl Rupp 81807171c5SHong Zhang Rt_dense->assembled = PETSC_TRUE; 82807171c5SHong Zhang rart->Rt = Rt_dense; 83807171c5SHong Zhang 84807171c5SHong Zhang /* Create RARt_dense = R*A*Rt_dense */ 85807171c5SHong Zhang ierr = MatCreate(PETSC_COMM_SELF,&RARt_dense);CHKERRQ(ierr); 86807171c5SHong Zhang ierr = MatSetSizes(RARt_dense,(*C)->rmap->n,matcoloring->ncolors,(*C)->rmap->n,matcoloring->ncolors);CHKERRQ(ierr); 87807171c5SHong Zhang ierr = MatSetType(RARt_dense,MATSEQDENSE);CHKERRQ(ierr); 880298fd71SBarry Smith ierr = MatSeqDenseSetPreallocation(RARt_dense,NULL);CHKERRQ(ierr); 892205254eSKarl Rupp 90807171c5SHong Zhang rart->RARt = RARt_dense; 91807171c5SHong Zhang 92807171c5SHong Zhang /* Allocate work array to store columns of A*R^T used in MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense() */ 93785e854fSJed Brown ierr = PetscMalloc1(A->rmap->n*4,&rart->work);CHKERRQ(ierr); 94807171c5SHong Zhang 95807171c5SHong Zhang rart->destroy = (*C)->ops->destroy; 96807171c5SHong Zhang (*C)->ops->destroy = MatDestroy_SeqAIJ_RARt; 97807171c5SHong Zhang 98807171c5SHong Zhang /* clean up */ 99807171c5SHong Zhang ierr = MatRestoreSymbolicTranspose_SeqAIJ(R,&rti,&rtj);CHKERRQ(ierr); 100807171c5SHong Zhang ierr = MatDestroy(&P);CHKERRQ(ierr); 101807171c5SHong Zhang 102807171c5SHong Zhang #if defined(PETSC_USE_INFO) 1031ce71dffSSatish Balay { 104807171c5SHong Zhang PetscReal density= (PetscReal)(c->nz)/(RARt_dense->rmap->n*RARt_dense->cmap->n); 1052cfd4408SHong Zhang ierr = PetscInfo(*C,"C=R*(A*Rt) via coloring C - use sparse-dense inner products\n");CHKERRQ(ierr); 1062cfd4408SHong Zhang ierr = PetscInfo6(*C,"RARt_den %D %D; Rt %D %D (RARt->nz %D)/(m*ncolors)=%g\n",RARt_dense->rmap->n,RARt_dense->cmap->n,R->cmap->n,R->rmap->n,c->nz,density);CHKERRQ(ierr); 1071ce71dffSSatish Balay } 108807171c5SHong Zhang #endif 109807171c5SHong Zhang PetscFunctionReturn(0); 110807171c5SHong Zhang } 111807171c5SHong Zhang 112807171c5SHong Zhang /* 113807171c5SHong Zhang RAB = R * A * B, R and A in seqaij format, B in dense format; 114807171c5SHong Zhang */ 115807171c5SHong Zhang PetscErrorCode MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense(Mat R,Mat A,Mat B,Mat RAB,PetscScalar *work) 116807171c5SHong Zhang { 117807171c5SHong Zhang Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*r=(Mat_SeqAIJ*)R->data; 118807171c5SHong Zhang PetscErrorCode ierr; 119807171c5SHong Zhang PetscScalar *b,r1,r2,r3,r4,*b1,*b2,*b3,*b4; 120807171c5SHong Zhang MatScalar *aa,*ra; 121807171c5SHong Zhang PetscInt cn =B->cmap->n,bm=B->rmap->n,col,i,j,n,*ai=a->i,*aj,am=A->rmap->n; 122807171c5SHong Zhang PetscInt am2=2*am,am3=3*am,bm4=4*bm; 123807171c5SHong Zhang PetscScalar *d,*c,*c2,*c3,*c4; 124807171c5SHong Zhang PetscInt *rj,rm=R->rmap->n,dm=RAB->rmap->n,dn=RAB->cmap->n; 125807171c5SHong Zhang PetscInt rm2=2*rm,rm3=3*rm,colrm; 126807171c5SHong Zhang 127807171c5SHong Zhang PetscFunctionBegin; 128807171c5SHong Zhang if (!dm || !dn) PetscFunctionReturn(0); 129807171c5SHong 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); 130807171c5SHong 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); 131807171c5SHong 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); 132807171c5SHong 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); 133807171c5SHong Zhang 134274010c0SHong Zhang { /* 135274010c0SHong Zhang This approach is not as good as original ones (will be removed later), but it reveals that 136274010c0SHong Zhang AB_den=A*B takes almost all execution time in R*A*B for src/ksp/ksp/examples/tutorials/ex56.c 137274010c0SHong Zhang */ 138274010c0SHong Zhang PetscBool via_matmatmult=PETSC_FALSE; 139c5929fdfSBarry Smith ierr = PetscOptionsGetBool(NULL,NULL,"-matrart_via_matmatmult",&via_matmatmult,NULL);CHKERRQ(ierr); 140274010c0SHong Zhang if (via_matmatmult) { 141274010c0SHong Zhang Mat AB_den; 142274010c0SHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqDense(A,B,0.0,&AB_den);CHKERRQ(ierr); 143274010c0SHong Zhang ierr = MatMatMultNumeric_SeqAIJ_SeqDense(A,B,AB_den);CHKERRQ(ierr); 144274010c0SHong Zhang ierr = MatMatMultNumeric_SeqAIJ_SeqDense(R,AB_den,RAB);CHKERRQ(ierr); 145274010c0SHong Zhang ierr = MatDestroy(&AB_den);CHKERRQ(ierr); 146274010c0SHong Zhang PetscFunctionReturn(0); 147274010c0SHong Zhang } 148274010c0SHong Zhang } 149274010c0SHong Zhang 1508c778c55SBarry Smith ierr = MatDenseGetArray(B,&b);CHKERRQ(ierr); 1518c778c55SBarry Smith ierr = MatDenseGetArray(RAB,&d);CHKERRQ(ierr); 152807171c5SHong Zhang b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm; 153807171c5SHong Zhang c = work; c2 = c + am; c3 = c2 + am; c4 = c3 + am; 154807171c5SHong Zhang for (col=0; col<cn-4; col += 4) { /* over columns of C */ 155807171c5SHong Zhang for (i=0; i<am; i++) { /* over rows of A in those columns */ 156807171c5SHong Zhang r1 = r2 = r3 = r4 = 0.0; 157807171c5SHong Zhang n = ai[i+1] - ai[i]; 158807171c5SHong Zhang aj = a->j + ai[i]; 159807171c5SHong Zhang aa = a->a + ai[i]; 160807171c5SHong Zhang for (j=0; j<n; j++) { 161807171c5SHong Zhang r1 += (*aa)*b1[*aj]; 162807171c5SHong Zhang r2 += (*aa)*b2[*aj]; 163807171c5SHong Zhang r3 += (*aa)*b3[*aj]; 164807171c5SHong Zhang r4 += (*aa++)*b4[*aj++]; 165807171c5SHong Zhang } 166807171c5SHong Zhang c[i] = r1; 167807171c5SHong Zhang c[am + i] = r2; 168807171c5SHong Zhang c[am2 + i] = r3; 169807171c5SHong Zhang c[am3 + i] = r4; 170807171c5SHong Zhang } 171807171c5SHong Zhang b1 += bm4; 172807171c5SHong Zhang b2 += bm4; 173807171c5SHong Zhang b3 += bm4; 174807171c5SHong Zhang b4 += bm4; 175807171c5SHong Zhang 176807171c5SHong Zhang /* RAB[:,col] = R*C[:,col] */ 177807171c5SHong Zhang colrm = col*rm; 178807171c5SHong Zhang for (i=0; i<rm; i++) { /* over rows of R in those columns */ 179807171c5SHong Zhang r1 = r2 = r3 = r4 = 0.0; 180807171c5SHong Zhang n = r->i[i+1] - r->i[i]; 181807171c5SHong Zhang rj = r->j + r->i[i]; 182807171c5SHong Zhang ra = r->a + r->i[i]; 183807171c5SHong Zhang for (j=0; j<n; j++) { 184807171c5SHong Zhang r1 += (*ra)*c[*rj]; 185807171c5SHong Zhang r2 += (*ra)*c2[*rj]; 186807171c5SHong Zhang r3 += (*ra)*c3[*rj]; 187807171c5SHong Zhang r4 += (*ra++)*c4[*rj++]; 188807171c5SHong Zhang } 189807171c5SHong Zhang d[colrm + i] = r1; 190807171c5SHong Zhang d[colrm + rm + i] = r2; 191807171c5SHong Zhang d[colrm + rm2 + i] = r3; 192807171c5SHong Zhang d[colrm + rm3 + i] = r4; 193807171c5SHong Zhang } 194807171c5SHong Zhang } 195807171c5SHong Zhang for (; col<cn; col++) { /* over extra columns of C */ 196807171c5SHong Zhang for (i=0; i<am; i++) { /* over rows of A in those columns */ 197807171c5SHong Zhang r1 = 0.0; 198807171c5SHong Zhang n = a->i[i+1] - a->i[i]; 199807171c5SHong Zhang aj = a->j + a->i[i]; 200807171c5SHong Zhang aa = a->a + a->i[i]; 201807171c5SHong Zhang for (j=0; j<n; j++) { 202807171c5SHong Zhang r1 += (*aa++)*b1[*aj++]; 203807171c5SHong Zhang } 204807171c5SHong Zhang c[i] = r1; 205807171c5SHong Zhang } 206807171c5SHong Zhang b1 += bm; 207807171c5SHong Zhang 208807171c5SHong Zhang for (i=0; i<rm; i++) { /* over rows of R in those columns */ 209807171c5SHong Zhang r1 = 0.0; 210807171c5SHong Zhang n = r->i[i+1] - r->i[i]; 211807171c5SHong Zhang rj = r->j + r->i[i]; 212807171c5SHong Zhang ra = r->a + r->i[i]; 213807171c5SHong Zhang for (j=0; j<n; j++) { 214807171c5SHong Zhang r1 += (*ra++)*c[*rj++]; 215807171c5SHong Zhang } 216807171c5SHong Zhang d[col*rm + i] = r1; 217807171c5SHong Zhang } 218807171c5SHong Zhang } 219807171c5SHong Zhang ierr = PetscLogFlops(cn*2.0*(a->nz + r->nz));CHKERRQ(ierr); 220807171c5SHong Zhang 2218c778c55SBarry Smith ierr = MatDenseRestoreArray(B,&b);CHKERRQ(ierr); 2228c778c55SBarry Smith ierr = MatDenseRestoreArray(RAB,&d);CHKERRQ(ierr); 223807171c5SHong Zhang ierr = MatAssemblyBegin(RAB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 224807171c5SHong Zhang ierr = MatAssemblyEnd(RAB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 225807171c5SHong Zhang PetscFunctionReturn(0); 226807171c5SHong Zhang } 227807171c5SHong Zhang 22855bea0ebSHong Zhang PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ_colorrart(Mat A,Mat R,Mat C) 229807171c5SHong Zhang { 230807171c5SHong Zhang PetscErrorCode ierr; 231257c235dSHong Zhang Mat_SeqAIJ *c = (Mat_SeqAIJ*)C->data; 232257c235dSHong Zhang Mat_RARt *rart=c->rart; 233807171c5SHong Zhang MatTransposeColoring matcoloring; 2348d0a38e4SHong Zhang Mat Rt,RARt; 235807171c5SHong Zhang 236807171c5SHong Zhang PetscFunctionBegin; 237807171c5SHong Zhang /* Get dense Rt by Apply MatTransposeColoring to R */ 238807171c5SHong Zhang matcoloring = rart->matcoloring; 239807171c5SHong Zhang Rt = rart->Rt; 240807171c5SHong Zhang ierr = MatTransColoringApplySpToDen(matcoloring,R,Rt);CHKERRQ(ierr); 241807171c5SHong Zhang 24250647e95SHong Zhang /* Get dense RARt = R*A*Rt -- dominates! */ 243807171c5SHong Zhang RARt = rart->RARt; 244807171c5SHong Zhang ierr = MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense(R,A,Rt,RARt,rart->work);CHKERRQ(ierr); 245807171c5SHong Zhang 246807171c5SHong Zhang /* Recover C from C_dense */ 247807171c5SHong Zhang ierr = MatTransColoringApplyDenToSp(matcoloring,RARt,C);CHKERRQ(ierr); 248807171c5SHong Zhang PetscFunctionReturn(0); 249807171c5SHong Zhang } 250807171c5SHong Zhang 25155bea0ebSHong Zhang PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ_matmattransposemult(Mat A,Mat R,PetscReal fill,Mat *C) 252807171c5SHong Zhang { 253807171c5SHong Zhang PetscErrorCode ierr; 25455bea0ebSHong Zhang Mat ARt,RARt; 2554fa85fdeSHong Zhang Mat_SeqAIJ *c; 2564fa85fdeSHong Zhang Mat_RARt *rart; 2574fa85fdeSHong Zhang 25895a72cc5SHong Zhang PetscFunctionBegin; 259b3b4fc5aSHong Zhang /* must use '-mat_no_inode' with '-matmattransmult_color 1' - do not knwo why? */ 26036104f73SHong Zhang ierr = MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(A,R,fill,&ARt);CHKERRQ(ierr); 26136104f73SHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(R,ARt,fill,&RARt);CHKERRQ(ierr); 2623b1b9624SHong Zhang *C = RARt; 26355bea0ebSHong Zhang RARt->ops->rartnumeric = MatRARtNumeric_SeqAIJ_SeqAIJ_matmattransposemult; 26455bea0ebSHong Zhang 265b00a9115SJed Brown ierr = PetscNew(&rart);CHKERRQ(ierr); 2663b1b9624SHong Zhang c = (Mat_SeqAIJ*)(*C)->data; 2673b1b9624SHong Zhang c->rart = rart; 2683b1b9624SHong Zhang rart->ARt = ARt; 2693b1b9624SHong Zhang rart->destroy = RARt->ops->destroy; 2703b1b9624SHong Zhang RARt->ops->destroy = MatDestroy_SeqAIJ_RARt; 2712cfd4408SHong Zhang #if defined(PETSC_USE_INFO) 2722cfd4408SHong 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); 2732cfd4408SHong Zhang #endif 27455bea0ebSHong Zhang PetscFunctionReturn(0); 27536104f73SHong Zhang } 27655bea0ebSHong Zhang 27755bea0ebSHong Zhang PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ_matmattransposemult(Mat A,Mat R,Mat C) 27855bea0ebSHong Zhang { 27955bea0ebSHong Zhang PetscErrorCode ierr; 28055bea0ebSHong Zhang Mat_SeqAIJ *c=(Mat_SeqAIJ*)C->data; 28155bea0ebSHong Zhang Mat_RARt *rart=c->rart; 28255bea0ebSHong Zhang Mat ARt=rart->ARt; 28355bea0ebSHong Zhang 28455bea0ebSHong Zhang PetscFunctionBegin; 285b3b4fc5aSHong Zhang ierr = MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(A,R,ARt);CHKERRQ(ierr); /* dominate! */ 28655bea0ebSHong Zhang ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ(R,ARt,C);CHKERRQ(ierr); 28755bea0ebSHong Zhang PetscFunctionReturn(0); 288807171c5SHong Zhang } 28955bea0ebSHong Zhang 29055bea0ebSHong Zhang PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat R,PetscReal fill,Mat *C) 29155bea0ebSHong Zhang { 29255bea0ebSHong Zhang PetscErrorCode ierr; 29395a72cc5SHong Zhang Mat Rt; 29455bea0ebSHong Zhang Mat_SeqAIJ *c; 29555bea0ebSHong Zhang Mat_RARt *rart; 29655bea0ebSHong Zhang 29755bea0ebSHong Zhang PetscFunctionBegin; 2985008f5a7SHong Zhang ierr = MatTranspose_SeqAIJ(R,MAT_INITIAL_MATRIX,&Rt);CHKERRQ(ierr); 2995008f5a7SHong Zhang ierr = MatMatMatMultSymbolic_SeqAIJ_SeqAIJ_SeqAIJ(R,A,Rt,fill,C);CHKERRQ(ierr); 30095a72cc5SHong Zhang 301b00a9115SJed Brown ierr = PetscNew(&rart);CHKERRQ(ierr); 30295a72cc5SHong Zhang rart->Rt = Rt; 30395a72cc5SHong Zhang c = (Mat_SeqAIJ*)(*C)->data; 30495a72cc5SHong Zhang c->rart = rart; 30595a72cc5SHong Zhang rart->destroy = (*C)->ops->destroy; 30695a72cc5SHong Zhang (*C)->ops->destroy = MatDestroy_SeqAIJ_RARt; 30755bea0ebSHong Zhang (*C)->ops->rartnumeric = MatRARtNumeric_SeqAIJ_SeqAIJ; 3082cfd4408SHong Zhang #if defined(PETSC_USE_INFO) 3092cfd4408SHong Zhang ierr = PetscInfo(*C,"Use Rt=R^T and C=R*A*Rt via MatMatMatMult() to avoid sparse inner products\n");CHKERRQ(ierr); 3102cfd4408SHong Zhang #endif 31155bea0ebSHong Zhang PetscFunctionReturn(0); 31255bea0ebSHong Zhang } 31355bea0ebSHong Zhang 31455bea0ebSHong Zhang PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ(Mat A,Mat R,Mat C) 31555bea0ebSHong Zhang { 31655bea0ebSHong Zhang PetscErrorCode ierr; 31755bea0ebSHong Zhang Mat_SeqAIJ *c = (Mat_SeqAIJ*)C->data; 31855bea0ebSHong Zhang Mat_RARt *rart = c->rart; 31955bea0ebSHong Zhang Mat Rt = rart->Rt; 32055bea0ebSHong Zhang 32155bea0ebSHong Zhang PetscFunctionBegin; 32255bea0ebSHong Zhang ierr = MatTranspose_SeqAIJ(R,MAT_REUSE_MATRIX,&Rt);CHKERRQ(ierr); 32355bea0ebSHong Zhang ierr = MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqAIJ(R,A,Rt,C);CHKERRQ(ierr); 32455bea0ebSHong Zhang PetscFunctionReturn(0); 32555bea0ebSHong Zhang } 32655bea0ebSHong Zhang 32755bea0ebSHong Zhang PetscErrorCode MatRARt_SeqAIJ_SeqAIJ(Mat A,Mat R,MatReuse scall,PetscReal fill,Mat *C) 32855bea0ebSHong Zhang { 32955bea0ebSHong Zhang PetscErrorCode ierr; 33055bea0ebSHong Zhang const char *algTypes[3] = {"matmatmatmult","matmattransposemult","coloring_rart"}; 33155bea0ebSHong Zhang PetscInt alg=0; /* set default algorithm */ 33255bea0ebSHong Zhang 33355bea0ebSHong Zhang PetscFunctionBegin; 33455bea0ebSHong Zhang if (scall == MAT_INITIAL_MATRIX) { 335*715a5346SHong Zhang ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"MatRARt","Mat");CHKERRQ(ierr); 33655bea0ebSHong Zhang ierr = PetscOptionsEList("-matrart_via","Algorithmic approach","MatRARt",algTypes,3,algTypes[0],&alg,NULL);CHKERRQ(ierr); 337b56132cbSHong Zhang ierr = PetscOptionsEnd();CHKERRQ(ierr); 338b56132cbSHong Zhang 33955bea0ebSHong Zhang ierr = PetscLogEventBegin(MAT_RARtSymbolic,A,R,0,0);CHKERRQ(ierr); 34055bea0ebSHong Zhang switch (alg) { 34155bea0ebSHong Zhang case 1: 34255bea0ebSHong Zhang /* via matmattransposemult: ARt=A*R^T, C=R*ARt - matrix coloring can be applied to A*R^T */ 34355bea0ebSHong Zhang ierr = MatRARtSymbolic_SeqAIJ_SeqAIJ_matmattransposemult(A,R,fill,C);CHKERRQ(ierr); 34455bea0ebSHong Zhang break; 34555bea0ebSHong Zhang case 2: 34655bea0ebSHong Zhang /* via coloring_rart: apply coloring C = R*A*R^T */ 34755bea0ebSHong Zhang ierr = MatRARtSymbolic_SeqAIJ_SeqAIJ_colorrart(A,R,fill,C);CHKERRQ(ierr); 34855bea0ebSHong Zhang break; 34955bea0ebSHong Zhang default: 35055bea0ebSHong Zhang /* via matmatmatmult: Rt=R^T, C=R*A*Rt - avoid inefficient sparse inner products */ 35155bea0ebSHong Zhang ierr = MatRARtSymbolic_SeqAIJ_SeqAIJ(A,R,fill,C);CHKERRQ(ierr); 35255bea0ebSHong Zhang break; 35355bea0ebSHong Zhang } 3545008f5a7SHong Zhang ierr = PetscLogEventEnd(MAT_RARtSymbolic,A,R,0,0);CHKERRQ(ierr); 3555008f5a7SHong Zhang } 3565008f5a7SHong Zhang 3575008f5a7SHong Zhang ierr = PetscLogEventBegin(MAT_RARtNumeric,A,R,0,0);CHKERRQ(ierr); 35855bea0ebSHong Zhang ierr = (*(*C)->ops->rartnumeric)(A,R,*C);CHKERRQ(ierr); 3595008f5a7SHong Zhang ierr = PetscLogEventEnd(MAT_RARtNumeric,A,R,0,0);CHKERRQ(ierr); 360807171c5SHong Zhang PetscFunctionReturn(0); 361807171c5SHong Zhang } 362