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 11807171c5SHong Zhang #undef __FUNCT__ 12257c235dSHong Zhang #define __FUNCT__ "MatDestroy_SeqAIJ_RARt" 13257c235dSHong Zhang PetscErrorCode MatDestroy_SeqAIJ_RARt(Mat A) 14807171c5SHong Zhang { 15807171c5SHong Zhang PetscErrorCode ierr; 16257c235dSHong Zhang Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 17257c235dSHong Zhang Mat_RARt *rart = a->rart; 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); 233b1b9624SHong Zhang ierr = MatDestroy(&rart->ARt);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__ 3555bea0ebSHong Zhang #define __FUNCT__ "MatRARtSymbolic_SeqAIJ_SeqAIJ_colorrart" 3655bea0ebSHong Zhang PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ_colorrart(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; 42335efc43SPeter Brune MatColoring coloring; 43807171c5SHong Zhang MatTransposeColoring matcoloring; 44807171c5SHong Zhang ISColoring iscoloring; 458d0a38e4SHong Zhang Mat Rt_dense,RARt_dense; 46807171c5SHong Zhang Mat_SeqAIJ *c; 47807171c5SHong Zhang 48807171c5SHong Zhang PetscFunctionBegin; 49807171c5SHong Zhang /* create symbolic P=Rt */ 50807171c5SHong Zhang ierr = MatGetSymbolicTranspose_SeqAIJ(R,&rti,&rtj);CHKERRQ(ierr); 510298fd71SBarry Smith ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,R->cmap->n,R->rmap->n,rti,rtj,NULL,&P);CHKERRQ(ierr); 52807171c5SHong Zhang 53807171c5SHong Zhang /* get symbolic C=Pt*A*P */ 544a1b09b7SHong Zhang ierr = MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy(A,P,fill,C);CHKERRQ(ierr); 55fd8e3dafSHong Zhang (*C)->rmap->bs = R->rmap->bs; 56fd8e3dafSHong Zhang (*C)->cmap->bs = R->rmap->bs; 5755bea0ebSHong Zhang (*C)->ops->rartnumeric = MatRARtNumeric_SeqAIJ_SeqAIJ_colorrart; 58807171c5SHong Zhang 59807171c5SHong Zhang /* create a supporting struct */ 60*b00a9115SJed Brown ierr = PetscNew(&rart);CHKERRQ(ierr); 61257c235dSHong Zhang c = (Mat_SeqAIJ*)(*C)->data; 62257c235dSHong Zhang c->rart = rart; 63807171c5SHong Zhang 6422e94b5dSHong Zhang /* ------ Use coloring ---------- */ 654d478ae7SHong Zhang /* inode causes memory problem, don't know why */ 664d478ae7SHong Zhang if (c->inode.use) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MAT_USE_INODES is not supported. Use '-mat_no_inode'"); 674d478ae7SHong Zhang 68807171c5SHong Zhang /* Create MatTransposeColoring from symbolic C=R*A*R^T */ 69335efc43SPeter Brune ierr = MatColoringCreate(*C,&coloring);CHKERRQ(ierr); 70335efc43SPeter Brune ierr = MatColoringSetDistance(coloring,2);CHKERRQ(ierr); 71335efc43SPeter Brune ierr = MatColoringSetType(coloring,MATCOLORINGSL);CHKERRQ(ierr); 72335efc43SPeter Brune ierr = MatColoringSetFromOptions(coloring);CHKERRQ(ierr); 73335efc43SPeter Brune ierr = MatColoringApply(coloring,&iscoloring);CHKERRQ(ierr); 74335efc43SPeter Brune ierr = MatColoringDestroy(&coloring);CHKERRQ(ierr); 75807171c5SHong Zhang ierr = MatTransposeColoringCreate(*C,iscoloring,&matcoloring);CHKERRQ(ierr); 762205254eSKarl Rupp 77807171c5SHong Zhang rart->matcoloring = matcoloring; 78807171c5SHong Zhang ierr = ISColoringDestroy(&iscoloring);CHKERRQ(ierr); 793b1b9624SHong Zhang 80807171c5SHong Zhang /* Create Rt_dense */ 81807171c5SHong Zhang ierr = MatCreate(PETSC_COMM_SELF,&Rt_dense);CHKERRQ(ierr); 82807171c5SHong Zhang ierr = MatSetSizes(Rt_dense,A->cmap->n,matcoloring->ncolors,A->cmap->n,matcoloring->ncolors);CHKERRQ(ierr); 83807171c5SHong Zhang ierr = MatSetType(Rt_dense,MATSEQDENSE);CHKERRQ(ierr); 840298fd71SBarry Smith ierr = MatSeqDenseSetPreallocation(Rt_dense,NULL);CHKERRQ(ierr); 852205254eSKarl Rupp 86807171c5SHong Zhang Rt_dense->assembled = PETSC_TRUE; 87807171c5SHong Zhang rart->Rt = Rt_dense; 88807171c5SHong Zhang 89807171c5SHong Zhang /* Create RARt_dense = R*A*Rt_dense */ 90807171c5SHong Zhang ierr = MatCreate(PETSC_COMM_SELF,&RARt_dense);CHKERRQ(ierr); 91807171c5SHong Zhang ierr = MatSetSizes(RARt_dense,(*C)->rmap->n,matcoloring->ncolors,(*C)->rmap->n,matcoloring->ncolors);CHKERRQ(ierr); 92807171c5SHong Zhang ierr = MatSetType(RARt_dense,MATSEQDENSE);CHKERRQ(ierr); 930298fd71SBarry Smith ierr = MatSeqDenseSetPreallocation(RARt_dense,NULL);CHKERRQ(ierr); 942205254eSKarl Rupp 95807171c5SHong Zhang rart->RARt = RARt_dense; 96807171c5SHong Zhang 97807171c5SHong Zhang /* Allocate work array to store columns of A*R^T used in MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense() */ 98785e854fSJed Brown ierr = PetscMalloc1(A->rmap->n*4,&rart->work);CHKERRQ(ierr); 99807171c5SHong Zhang 100807171c5SHong Zhang rart->destroy = (*C)->ops->destroy; 101807171c5SHong Zhang (*C)->ops->destroy = MatDestroy_SeqAIJ_RARt; 102807171c5SHong Zhang 103807171c5SHong Zhang /* clean up */ 104807171c5SHong Zhang ierr = MatRestoreSymbolicTranspose_SeqAIJ(R,&rti,&rtj);CHKERRQ(ierr); 105807171c5SHong Zhang ierr = MatDestroy(&P);CHKERRQ(ierr); 106807171c5SHong Zhang 107807171c5SHong Zhang #if defined(PETSC_USE_INFO) 1081ce71dffSSatish Balay { 109807171c5SHong Zhang PetscReal density= (PetscReal)(c->nz)/(RARt_dense->rmap->n*RARt_dense->cmap->n); 1102cfd4408SHong Zhang ierr = PetscInfo(*C,"C=R*(A*Rt) via coloring C - use sparse-dense inner products\n");CHKERRQ(ierr); 1112cfd4408SHong 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); 1121ce71dffSSatish Balay } 113807171c5SHong Zhang #endif 114807171c5SHong Zhang PetscFunctionReturn(0); 115807171c5SHong Zhang } 116807171c5SHong Zhang 117807171c5SHong Zhang /* 118807171c5SHong Zhang RAB = R * A * B, R and A in seqaij format, B in dense format; 119807171c5SHong Zhang */ 120807171c5SHong Zhang #undef __FUNCT__ 121807171c5SHong Zhang #define __FUNCT__ "MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense" 122807171c5SHong Zhang PetscErrorCode MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense(Mat R,Mat A,Mat B,Mat RAB,PetscScalar *work) 123807171c5SHong Zhang { 124807171c5SHong Zhang Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*r=(Mat_SeqAIJ*)R->data; 125807171c5SHong Zhang PetscErrorCode ierr; 126807171c5SHong Zhang PetscScalar *b,r1,r2,r3,r4,*b1,*b2,*b3,*b4; 127807171c5SHong Zhang MatScalar *aa,*ra; 128807171c5SHong Zhang PetscInt cn =B->cmap->n,bm=B->rmap->n,col,i,j,n,*ai=a->i,*aj,am=A->rmap->n; 129807171c5SHong Zhang PetscInt am2=2*am,am3=3*am,bm4=4*bm; 130807171c5SHong Zhang PetscScalar *d,*c,*c2,*c3,*c4; 131807171c5SHong Zhang PetscInt *rj,rm=R->rmap->n,dm=RAB->rmap->n,dn=RAB->cmap->n; 132807171c5SHong Zhang PetscInt rm2=2*rm,rm3=3*rm,colrm; 133807171c5SHong Zhang 134807171c5SHong Zhang PetscFunctionBegin; 135807171c5SHong Zhang if (!dm || !dn) PetscFunctionReturn(0); 136807171c5SHong 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); 137807171c5SHong 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); 138807171c5SHong 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); 139807171c5SHong 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); 140807171c5SHong Zhang 141274010c0SHong Zhang { /* 142274010c0SHong Zhang This approach is not as good as original ones (will be removed later), but it reveals that 143274010c0SHong Zhang AB_den=A*B takes almost all execution time in R*A*B for src/ksp/ksp/examples/tutorials/ex56.c 144274010c0SHong Zhang */ 145274010c0SHong Zhang PetscBool via_matmatmult=PETSC_FALSE; 146274010c0SHong Zhang ierr = PetscOptionsGetBool(NULL,"-matrart_via_matmatmult",&via_matmatmult,NULL);CHKERRQ(ierr); 147274010c0SHong Zhang if (via_matmatmult) { 148274010c0SHong Zhang Mat AB_den; 149274010c0SHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqDense(A,B,0.0,&AB_den);CHKERRQ(ierr); 150274010c0SHong Zhang ierr = MatMatMultNumeric_SeqAIJ_SeqDense(A,B,AB_den);CHKERRQ(ierr); 151274010c0SHong Zhang ierr = MatMatMultNumeric_SeqAIJ_SeqDense(R,AB_den,RAB);CHKERRQ(ierr); 152274010c0SHong Zhang ierr = MatDestroy(&AB_den);CHKERRQ(ierr); 153274010c0SHong Zhang PetscFunctionReturn(0); 154274010c0SHong Zhang } 155274010c0SHong Zhang } 156274010c0SHong Zhang 1578c778c55SBarry Smith ierr = MatDenseGetArray(B,&b);CHKERRQ(ierr); 1588c778c55SBarry Smith ierr = MatDenseGetArray(RAB,&d);CHKERRQ(ierr); 159807171c5SHong Zhang b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm; 160807171c5SHong Zhang c = work; c2 = c + am; c3 = c2 + am; c4 = c3 + am; 161807171c5SHong Zhang for (col=0; col<cn-4; col += 4) { /* over columns of C */ 162807171c5SHong Zhang for (i=0; i<am; i++) { /* over rows of A in those columns */ 163807171c5SHong Zhang r1 = r2 = r3 = r4 = 0.0; 164807171c5SHong Zhang n = ai[i+1] - ai[i]; 165807171c5SHong Zhang aj = a->j + ai[i]; 166807171c5SHong Zhang aa = a->a + ai[i]; 167807171c5SHong Zhang for (j=0; j<n; j++) { 168807171c5SHong Zhang r1 += (*aa)*b1[*aj]; 169807171c5SHong Zhang r2 += (*aa)*b2[*aj]; 170807171c5SHong Zhang r3 += (*aa)*b3[*aj]; 171807171c5SHong Zhang r4 += (*aa++)*b4[*aj++]; 172807171c5SHong Zhang } 173807171c5SHong Zhang c[i] = r1; 174807171c5SHong Zhang c[am + i] = r2; 175807171c5SHong Zhang c[am2 + i] = r3; 176807171c5SHong Zhang c[am3 + i] = r4; 177807171c5SHong Zhang } 178807171c5SHong Zhang b1 += bm4; 179807171c5SHong Zhang b2 += bm4; 180807171c5SHong Zhang b3 += bm4; 181807171c5SHong Zhang b4 += bm4; 182807171c5SHong Zhang 183807171c5SHong Zhang /* RAB[:,col] = R*C[:,col] */ 184807171c5SHong Zhang colrm = col*rm; 185807171c5SHong Zhang for (i=0; i<rm; i++) { /* over rows of R in those columns */ 186807171c5SHong Zhang r1 = r2 = r3 = r4 = 0.0; 187807171c5SHong Zhang n = r->i[i+1] - r->i[i]; 188807171c5SHong Zhang rj = r->j + r->i[i]; 189807171c5SHong Zhang ra = r->a + r->i[i]; 190807171c5SHong Zhang for (j=0; j<n; j++) { 191807171c5SHong Zhang r1 += (*ra)*c[*rj]; 192807171c5SHong Zhang r2 += (*ra)*c2[*rj]; 193807171c5SHong Zhang r3 += (*ra)*c3[*rj]; 194807171c5SHong Zhang r4 += (*ra++)*c4[*rj++]; 195807171c5SHong Zhang } 196807171c5SHong Zhang d[colrm + i] = r1; 197807171c5SHong Zhang d[colrm + rm + i] = r2; 198807171c5SHong Zhang d[colrm + rm2 + i] = r3; 199807171c5SHong Zhang d[colrm + rm3 + i] = r4; 200807171c5SHong Zhang } 201807171c5SHong Zhang } 202807171c5SHong Zhang for (; col<cn; col++) { /* over extra columns of C */ 203807171c5SHong Zhang for (i=0; i<am; i++) { /* over rows of A in those columns */ 204807171c5SHong Zhang r1 = 0.0; 205807171c5SHong Zhang n = a->i[i+1] - a->i[i]; 206807171c5SHong Zhang aj = a->j + a->i[i]; 207807171c5SHong Zhang aa = a->a + a->i[i]; 208807171c5SHong Zhang for (j=0; j<n; j++) { 209807171c5SHong Zhang r1 += (*aa++)*b1[*aj++]; 210807171c5SHong Zhang } 211807171c5SHong Zhang c[i] = r1; 212807171c5SHong Zhang } 213807171c5SHong Zhang b1 += bm; 214807171c5SHong Zhang 215807171c5SHong Zhang for (i=0; i<rm; i++) { /* over rows of R in those columns */ 216807171c5SHong Zhang r1 = 0.0; 217807171c5SHong Zhang n = r->i[i+1] - r->i[i]; 218807171c5SHong Zhang rj = r->j + r->i[i]; 219807171c5SHong Zhang ra = r->a + r->i[i]; 220807171c5SHong Zhang for (j=0; j<n; j++) { 221807171c5SHong Zhang r1 += (*ra++)*c[*rj++]; 222807171c5SHong Zhang } 223807171c5SHong Zhang d[col*rm + i] = r1; 224807171c5SHong Zhang } 225807171c5SHong Zhang } 226807171c5SHong Zhang ierr = PetscLogFlops(cn*2.0*(a->nz + r->nz));CHKERRQ(ierr); 227807171c5SHong Zhang 2288c778c55SBarry Smith ierr = MatDenseRestoreArray(B,&b);CHKERRQ(ierr); 2298c778c55SBarry Smith ierr = MatDenseRestoreArray(RAB,&d);CHKERRQ(ierr); 230807171c5SHong Zhang ierr = MatAssemblyBegin(RAB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 231807171c5SHong Zhang ierr = MatAssemblyEnd(RAB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 232807171c5SHong Zhang PetscFunctionReturn(0); 233807171c5SHong Zhang } 234807171c5SHong Zhang 235807171c5SHong Zhang #undef __FUNCT__ 23655bea0ebSHong Zhang #define __FUNCT__ "MatRARtNumeric_SeqAIJ_SeqAIJ_colorrart" 23755bea0ebSHong Zhang PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ_colorrart(Mat A,Mat R,Mat C) 238807171c5SHong Zhang { 239807171c5SHong Zhang PetscErrorCode ierr; 240257c235dSHong Zhang Mat_SeqAIJ *c = (Mat_SeqAIJ*)C->data; 241257c235dSHong Zhang Mat_RARt *rart=c->rart; 242807171c5SHong Zhang MatTransposeColoring matcoloring; 2438d0a38e4SHong Zhang Mat Rt,RARt; 244807171c5SHong Zhang 245807171c5SHong Zhang PetscFunctionBegin; 246807171c5SHong Zhang /* Get dense Rt by Apply MatTransposeColoring to R */ 247807171c5SHong Zhang matcoloring = rart->matcoloring; 248807171c5SHong Zhang Rt = rart->Rt; 249807171c5SHong Zhang ierr = MatTransColoringApplySpToDen(matcoloring,R,Rt);CHKERRQ(ierr); 250807171c5SHong Zhang 25150647e95SHong Zhang /* Get dense RARt = R*A*Rt -- dominates! */ 252807171c5SHong Zhang RARt = rart->RARt; 253807171c5SHong Zhang ierr = MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense(R,A,Rt,RARt,rart->work);CHKERRQ(ierr); 254807171c5SHong Zhang 255807171c5SHong Zhang /* Recover C from C_dense */ 256807171c5SHong Zhang ierr = MatTransColoringApplyDenToSp(matcoloring,RARt,C);CHKERRQ(ierr); 257807171c5SHong Zhang PetscFunctionReturn(0); 258807171c5SHong Zhang } 259807171c5SHong Zhang 260807171c5SHong Zhang #undef __FUNCT__ 26155bea0ebSHong Zhang #define __FUNCT__ "MatRARtSymbolic_SeqAIJ_SeqAIJ_matmattransposemult" 26255bea0ebSHong Zhang PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ_matmattransposemult(Mat A,Mat R,PetscReal fill,Mat *C) 263807171c5SHong Zhang { 264807171c5SHong Zhang PetscErrorCode ierr; 26555bea0ebSHong Zhang Mat ARt,RARt; 2664fa85fdeSHong Zhang Mat_SeqAIJ *c; 2674fa85fdeSHong Zhang Mat_RARt *rart; 2684fa85fdeSHong Zhang 26995a72cc5SHong Zhang PetscFunctionBegin; 270b3b4fc5aSHong Zhang /* must use '-mat_no_inode' with '-matmattransmult_color 1' - do not knwo why? */ 27136104f73SHong Zhang ierr = MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(A,R,fill,&ARt);CHKERRQ(ierr); 27236104f73SHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(R,ARt,fill,&RARt);CHKERRQ(ierr); 2733b1b9624SHong Zhang *C = RARt; 27455bea0ebSHong Zhang RARt->ops->rartnumeric = MatRARtNumeric_SeqAIJ_SeqAIJ_matmattransposemult; 27555bea0ebSHong Zhang 276*b00a9115SJed Brown ierr = PetscNew(&rart);CHKERRQ(ierr); 2773b1b9624SHong Zhang c = (Mat_SeqAIJ*)(*C)->data; 2783b1b9624SHong Zhang c->rart = rart; 2793b1b9624SHong Zhang rart->ARt = ARt; 2803b1b9624SHong Zhang rart->destroy = RARt->ops->destroy; 2813b1b9624SHong Zhang RARt->ops->destroy = MatDestroy_SeqAIJ_RARt; 2822cfd4408SHong Zhang #if defined(PETSC_USE_INFO) 2832cfd4408SHong 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); 2842cfd4408SHong Zhang #endif 28555bea0ebSHong Zhang PetscFunctionReturn(0); 28636104f73SHong Zhang } 28755bea0ebSHong Zhang 28855bea0ebSHong Zhang #undef __FUNCT__ 28955bea0ebSHong Zhang #define __FUNCT__ "MatRARtNumeric_SeqAIJ_SeqAIJ_matmattransposemult" 29055bea0ebSHong Zhang PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ_matmattransposemult(Mat A,Mat R,Mat C) 29155bea0ebSHong Zhang { 29255bea0ebSHong Zhang PetscErrorCode ierr; 29355bea0ebSHong Zhang Mat_SeqAIJ *c=(Mat_SeqAIJ*)C->data; 29455bea0ebSHong Zhang Mat_RARt *rart=c->rart; 29555bea0ebSHong Zhang Mat ARt=rart->ARt; 29655bea0ebSHong Zhang 29755bea0ebSHong Zhang PetscFunctionBegin; 298b3b4fc5aSHong Zhang ierr = MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(A,R,ARt);CHKERRQ(ierr); /* dominate! */ 29955bea0ebSHong Zhang ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ(R,ARt,C);CHKERRQ(ierr); 30055bea0ebSHong Zhang PetscFunctionReturn(0); 301807171c5SHong Zhang } 30255bea0ebSHong Zhang 30355bea0ebSHong Zhang #undef __FUNCT__ 30455bea0ebSHong Zhang #define __FUNCT__ "MatRARtSymbolic_SeqAIJ_SeqAIJ" 30555bea0ebSHong Zhang PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat R,PetscReal fill,Mat *C) 30655bea0ebSHong Zhang { 30755bea0ebSHong Zhang PetscErrorCode ierr; 30895a72cc5SHong Zhang Mat Rt; 30955bea0ebSHong Zhang Mat_SeqAIJ *c; 31055bea0ebSHong Zhang Mat_RARt *rart; 31155bea0ebSHong Zhang 31255bea0ebSHong Zhang PetscFunctionBegin; 3135008f5a7SHong Zhang ierr = MatTranspose_SeqAIJ(R,MAT_INITIAL_MATRIX,&Rt);CHKERRQ(ierr); 3145008f5a7SHong Zhang ierr = MatMatMatMultSymbolic_SeqAIJ_SeqAIJ_SeqAIJ(R,A,Rt,fill,C);CHKERRQ(ierr); 31595a72cc5SHong Zhang 316*b00a9115SJed Brown ierr = PetscNew(&rart);CHKERRQ(ierr); 31795a72cc5SHong Zhang rart->Rt = Rt; 31895a72cc5SHong Zhang c = (Mat_SeqAIJ*)(*C)->data; 31995a72cc5SHong Zhang c->rart = rart; 32095a72cc5SHong Zhang rart->destroy = (*C)->ops->destroy; 32195a72cc5SHong Zhang (*C)->ops->destroy = MatDestroy_SeqAIJ_RARt; 32255bea0ebSHong Zhang (*C)->ops->rartnumeric = MatRARtNumeric_SeqAIJ_SeqAIJ; 3232cfd4408SHong Zhang #if defined(PETSC_USE_INFO) 3242cfd4408SHong Zhang ierr = PetscInfo(*C,"Use Rt=R^T and C=R*A*Rt via MatMatMatMult() to avoid sparse inner products\n");CHKERRQ(ierr); 3252cfd4408SHong Zhang #endif 32655bea0ebSHong Zhang PetscFunctionReturn(0); 32755bea0ebSHong Zhang } 32855bea0ebSHong Zhang 32955bea0ebSHong Zhang #undef __FUNCT__ 33055bea0ebSHong Zhang #define __FUNCT__ "MatRARtNumeric_SeqAIJ_SeqAIJ" 33155bea0ebSHong Zhang PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ(Mat A,Mat R,Mat C) 33255bea0ebSHong Zhang { 33355bea0ebSHong Zhang PetscErrorCode ierr; 33455bea0ebSHong Zhang Mat_SeqAIJ *c = (Mat_SeqAIJ*)C->data; 33555bea0ebSHong Zhang Mat_RARt *rart = c->rart; 33655bea0ebSHong Zhang Mat Rt = rart->Rt; 33755bea0ebSHong Zhang 33855bea0ebSHong Zhang PetscFunctionBegin; 33955bea0ebSHong Zhang ierr = MatTranspose_SeqAIJ(R,MAT_REUSE_MATRIX,&Rt);CHKERRQ(ierr); 34055bea0ebSHong Zhang ierr = MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqAIJ(R,A,Rt,C);CHKERRQ(ierr); 34155bea0ebSHong Zhang PetscFunctionReturn(0); 34255bea0ebSHong Zhang } 34355bea0ebSHong Zhang 34455bea0ebSHong Zhang #undef __FUNCT__ 34555bea0ebSHong Zhang #define __FUNCT__ "MatRARt_SeqAIJ_SeqAIJ" 34655bea0ebSHong Zhang PetscErrorCode MatRARt_SeqAIJ_SeqAIJ(Mat A,Mat R,MatReuse scall,PetscReal fill,Mat *C) 34755bea0ebSHong Zhang { 34855bea0ebSHong Zhang PetscErrorCode ierr; 34955bea0ebSHong Zhang const char *algTypes[3] = {"matmatmatmult","matmattransposemult","coloring_rart"}; 35055bea0ebSHong Zhang PetscInt alg=0; /* set default algorithm */ 35155bea0ebSHong Zhang 35255bea0ebSHong Zhang PetscFunctionBegin; 35355bea0ebSHong Zhang if (scall == MAT_INITIAL_MATRIX) { 354b56132cbSHong Zhang ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr); 35555bea0ebSHong Zhang ierr = PetscOptionsEList("-matrart_via","Algorithmic approach","MatRARt",algTypes,3,algTypes[0],&alg,NULL);CHKERRQ(ierr); 356b56132cbSHong Zhang ierr = PetscOptionsEnd();CHKERRQ(ierr); 357b56132cbSHong Zhang 35855bea0ebSHong Zhang ierr = PetscLogEventBegin(MAT_RARtSymbolic,A,R,0,0);CHKERRQ(ierr); 35955bea0ebSHong Zhang switch (alg) { 36055bea0ebSHong Zhang case 1: 36155bea0ebSHong Zhang /* via matmattransposemult: ARt=A*R^T, C=R*ARt - matrix coloring can be applied to A*R^T */ 36255bea0ebSHong Zhang ierr = MatRARtSymbolic_SeqAIJ_SeqAIJ_matmattransposemult(A,R,fill,C);CHKERRQ(ierr); 36355bea0ebSHong Zhang break; 36455bea0ebSHong Zhang case 2: 36555bea0ebSHong Zhang /* via coloring_rart: apply coloring C = R*A*R^T */ 36655bea0ebSHong Zhang ierr = MatRARtSymbolic_SeqAIJ_SeqAIJ_colorrart(A,R,fill,C);CHKERRQ(ierr); 36755bea0ebSHong Zhang break; 36855bea0ebSHong Zhang default: 36955bea0ebSHong Zhang /* via matmatmatmult: Rt=R^T, C=R*A*Rt - avoid inefficient sparse inner products */ 37055bea0ebSHong Zhang ierr = MatRARtSymbolic_SeqAIJ_SeqAIJ(A,R,fill,C);CHKERRQ(ierr); 37155bea0ebSHong Zhang break; 37255bea0ebSHong Zhang } 3735008f5a7SHong Zhang ierr = PetscLogEventEnd(MAT_RARtSymbolic,A,R,0,0);CHKERRQ(ierr); 3745008f5a7SHong Zhang } 3755008f5a7SHong Zhang 3765008f5a7SHong Zhang ierr = PetscLogEventBegin(MAT_RARtNumeric,A,R,0,0);CHKERRQ(ierr); 37755bea0ebSHong Zhang ierr = (*(*C)->ops->rartnumeric)(A,R,*C);CHKERRQ(ierr); 3785008f5a7SHong Zhang ierr = PetscLogEventEnd(MAT_RARtNumeric,A,R,0,0);CHKERRQ(ierr); 379807171c5SHong Zhang PetscFunctionReturn(0); 380807171c5SHong Zhang } 381