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*/ 108563dfccSBarry Smith #include <petsctime.h> 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); 243b1b9624SHong Zhang ierr = MatDestroy(&rart->ARt);CHKERRQ(ierr); 25807171c5SHong Zhang ierr = PetscFree(rart->work);CHKERRQ(ierr); 26807171c5SHong Zhang 27807171c5SHong Zhang A->ops->destroy = rart->destroy; 28807171c5SHong Zhang if (A->ops->destroy) { 29807171c5SHong Zhang ierr = (*A->ops->destroy)(A);CHKERRQ(ierr); 30807171c5SHong Zhang } 31257c235dSHong Zhang ierr = PetscFree(rart);CHKERRQ(ierr); 32807171c5SHong Zhang PetscFunctionReturn(0); 33807171c5SHong Zhang } 34807171c5SHong Zhang 35807171c5SHong Zhang #undef __FUNCT__ 3655bea0ebSHong Zhang #define __FUNCT__ "MatRARtSymbolic_SeqAIJ_SeqAIJ_colorrart" 3755bea0ebSHong Zhang PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ_colorrart(Mat A,Mat R,PetscReal fill,Mat *C) 38807171c5SHong Zhang { 39807171c5SHong Zhang PetscErrorCode ierr; 40807171c5SHong Zhang Mat P; 41807171c5SHong Zhang PetscInt *rti,*rtj; 42807171c5SHong Zhang Mat_RARt *rart; 43*335efc43SPeter Brune MatColoring coloring; 44807171c5SHong Zhang MatTransposeColoring matcoloring; 45807171c5SHong Zhang ISColoring iscoloring; 468d0a38e4SHong Zhang Mat Rt_dense,RARt_dense; 47807171c5SHong Zhang Mat_SeqAIJ *c; 48807171c5SHong Zhang 49807171c5SHong Zhang PetscFunctionBegin; 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 */ 554a1b09b7SHong Zhang ierr = MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy(A,P,fill,C);CHKERRQ(ierr); 56fd8e3dafSHong Zhang (*C)->rmap->bs = R->rmap->bs; 57fd8e3dafSHong Zhang (*C)->cmap->bs = R->rmap->bs; 5855bea0ebSHong Zhang (*C)->ops->rartnumeric = MatRARtNumeric_SeqAIJ_SeqAIJ_colorrart; 59807171c5SHong Zhang 60807171c5SHong Zhang /* create a supporting struct */ 61807171c5SHong Zhang ierr = PetscNew(Mat_RARt,&rart);CHKERRQ(ierr); 62257c235dSHong Zhang c = (Mat_SeqAIJ*)(*C)->data; 63257c235dSHong Zhang c->rart = rart; 64807171c5SHong Zhang 6522e94b5dSHong Zhang /* ------ Use coloring ---------- */ 664d478ae7SHong Zhang /* inode causes memory problem, don't know why */ 674d478ae7SHong Zhang if (c->inode.use) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MAT_USE_INODES is not supported. Use '-mat_no_inode'"); 684d478ae7SHong Zhang 69807171c5SHong Zhang /* Create MatTransposeColoring from symbolic C=R*A*R^T */ 70*335efc43SPeter Brune ierr = MatColoringCreate(*C,&coloring);CHKERRQ(ierr); 71*335efc43SPeter Brune ierr = MatColoringSetDistance(coloring,2);CHKERRQ(ierr); 72*335efc43SPeter Brune ierr = MatColoringSetType(coloring,MATCOLORINGSL);CHKERRQ(ierr); 73*335efc43SPeter Brune ierr = MatColoringSetFromOptions(coloring);CHKERRQ(ierr); 74*335efc43SPeter Brune ierr = MatColoringApply(coloring,&iscoloring);CHKERRQ(ierr); 75*335efc43SPeter Brune ierr = MatColoringDestroy(&coloring);CHKERRQ(ierr); 76807171c5SHong Zhang ierr = MatTransposeColoringCreate(*C,iscoloring,&matcoloring);CHKERRQ(ierr); 772205254eSKarl Rupp 78807171c5SHong Zhang rart->matcoloring = matcoloring; 79807171c5SHong Zhang ierr = ISColoringDestroy(&iscoloring);CHKERRQ(ierr); 803b1b9624SHong Zhang 81807171c5SHong Zhang /* Create Rt_dense */ 82807171c5SHong Zhang ierr = MatCreate(PETSC_COMM_SELF,&Rt_dense);CHKERRQ(ierr); 83807171c5SHong Zhang ierr = MatSetSizes(Rt_dense,A->cmap->n,matcoloring->ncolors,A->cmap->n,matcoloring->ncolors);CHKERRQ(ierr); 84807171c5SHong Zhang ierr = MatSetType(Rt_dense,MATSEQDENSE);CHKERRQ(ierr); 850298fd71SBarry Smith ierr = MatSeqDenseSetPreallocation(Rt_dense,NULL);CHKERRQ(ierr); 862205254eSKarl Rupp 87807171c5SHong Zhang Rt_dense->assembled = PETSC_TRUE; 88807171c5SHong Zhang rart->Rt = Rt_dense; 89807171c5SHong Zhang 90807171c5SHong Zhang /* Create RARt_dense = R*A*Rt_dense */ 91807171c5SHong Zhang ierr = MatCreate(PETSC_COMM_SELF,&RARt_dense);CHKERRQ(ierr); 92807171c5SHong Zhang ierr = MatSetSizes(RARt_dense,(*C)->rmap->n,matcoloring->ncolors,(*C)->rmap->n,matcoloring->ncolors);CHKERRQ(ierr); 93807171c5SHong Zhang ierr = MatSetType(RARt_dense,MATSEQDENSE);CHKERRQ(ierr); 940298fd71SBarry Smith ierr = MatSeqDenseSetPreallocation(RARt_dense,NULL);CHKERRQ(ierr); 952205254eSKarl Rupp 96807171c5SHong Zhang rart->RARt = RARt_dense; 97807171c5SHong Zhang 98807171c5SHong Zhang /* Allocate work array to store columns of A*R^T used in MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense() */ 99807171c5SHong Zhang ierr = PetscMalloc(A->rmap->n*4*sizeof(PetscScalar),&rart->work);CHKERRQ(ierr); 100807171c5SHong Zhang 101807171c5SHong Zhang rart->destroy = (*C)->ops->destroy; 102807171c5SHong Zhang (*C)->ops->destroy = MatDestroy_SeqAIJ_RARt; 103807171c5SHong Zhang 104807171c5SHong Zhang /* clean up */ 105807171c5SHong Zhang ierr = MatRestoreSymbolicTranspose_SeqAIJ(R,&rti,&rtj);CHKERRQ(ierr); 106807171c5SHong Zhang ierr = MatDestroy(&P);CHKERRQ(ierr); 107807171c5SHong Zhang 108807171c5SHong Zhang #if defined(PETSC_USE_INFO) 1091ce71dffSSatish Balay { 110807171c5SHong Zhang PetscReal density= (PetscReal)(c->nz)/(RARt_dense->rmap->n*RARt_dense->cmap->n); 1112cfd4408SHong Zhang ierr = PetscInfo(*C,"C=R*(A*Rt) via coloring C - use sparse-dense inner products\n");CHKERRQ(ierr); 1122cfd4408SHong 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); 1131ce71dffSSatish Balay } 114807171c5SHong Zhang #endif 115807171c5SHong Zhang PetscFunctionReturn(0); 116807171c5SHong Zhang } 117807171c5SHong Zhang 118807171c5SHong Zhang /* 119807171c5SHong Zhang RAB = R * A * B, R and A in seqaij format, B in dense format; 120807171c5SHong Zhang */ 121807171c5SHong Zhang #undef __FUNCT__ 122807171c5SHong Zhang #define __FUNCT__ "MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense" 123807171c5SHong Zhang PetscErrorCode MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense(Mat R,Mat A,Mat B,Mat RAB,PetscScalar *work) 124807171c5SHong Zhang { 125807171c5SHong Zhang Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*r=(Mat_SeqAIJ*)R->data; 126807171c5SHong Zhang PetscErrorCode ierr; 127807171c5SHong Zhang PetscScalar *b,r1,r2,r3,r4,*b1,*b2,*b3,*b4; 128807171c5SHong Zhang MatScalar *aa,*ra; 129807171c5SHong Zhang PetscInt cn =B->cmap->n,bm=B->rmap->n,col,i,j,n,*ai=a->i,*aj,am=A->rmap->n; 130807171c5SHong Zhang PetscInt am2=2*am,am3=3*am,bm4=4*bm; 131807171c5SHong Zhang PetscScalar *d,*c,*c2,*c3,*c4; 132807171c5SHong Zhang PetscInt *rj,rm=R->rmap->n,dm=RAB->rmap->n,dn=RAB->cmap->n; 133807171c5SHong Zhang PetscInt rm2=2*rm,rm3=3*rm,colrm; 134807171c5SHong Zhang 135807171c5SHong Zhang PetscFunctionBegin; 136807171c5SHong Zhang if (!dm || !dn) PetscFunctionReturn(0); 137807171c5SHong 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); 138807171c5SHong 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); 139807171c5SHong 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); 140807171c5SHong 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); 141807171c5SHong Zhang 142274010c0SHong Zhang { /* 143274010c0SHong Zhang This approach is not as good as original ones (will be removed later), but it reveals that 144274010c0SHong Zhang AB_den=A*B takes almost all execution time in R*A*B for src/ksp/ksp/examples/tutorials/ex56.c 145274010c0SHong Zhang */ 146274010c0SHong Zhang PetscBool via_matmatmult=PETSC_FALSE; 147274010c0SHong Zhang ierr = PetscOptionsGetBool(NULL,"-matrart_via_matmatmult",&via_matmatmult,NULL);CHKERRQ(ierr); 148274010c0SHong Zhang if (via_matmatmult) { 149274010c0SHong Zhang Mat AB_den; 150274010c0SHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqDense(A,B,0.0,&AB_den);CHKERRQ(ierr); 151274010c0SHong Zhang ierr = MatMatMultNumeric_SeqAIJ_SeqDense(A,B,AB_den);CHKERRQ(ierr); 152274010c0SHong Zhang ierr = MatMatMultNumeric_SeqAIJ_SeqDense(R,AB_den,RAB);CHKERRQ(ierr); 153274010c0SHong Zhang ierr = MatDestroy(&AB_den);CHKERRQ(ierr); 154274010c0SHong Zhang PetscFunctionReturn(0); 155274010c0SHong Zhang } 156274010c0SHong Zhang } 157274010c0SHong Zhang 1588c778c55SBarry Smith ierr = MatDenseGetArray(B,&b);CHKERRQ(ierr); 1598c778c55SBarry Smith ierr = MatDenseGetArray(RAB,&d);CHKERRQ(ierr); 160807171c5SHong Zhang b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm; 161807171c5SHong Zhang c = work; c2 = c + am; c3 = c2 + am; c4 = c3 + am; 162807171c5SHong Zhang for (col=0; col<cn-4; col += 4) { /* over columns of C */ 163807171c5SHong Zhang for (i=0; i<am; i++) { /* over rows of A in those columns */ 164807171c5SHong Zhang r1 = r2 = r3 = r4 = 0.0; 165807171c5SHong Zhang n = ai[i+1] - ai[i]; 166807171c5SHong Zhang aj = a->j + ai[i]; 167807171c5SHong Zhang aa = a->a + ai[i]; 168807171c5SHong Zhang for (j=0; j<n; j++) { 169807171c5SHong Zhang r1 += (*aa)*b1[*aj]; 170807171c5SHong Zhang r2 += (*aa)*b2[*aj]; 171807171c5SHong Zhang r3 += (*aa)*b3[*aj]; 172807171c5SHong Zhang r4 += (*aa++)*b4[*aj++]; 173807171c5SHong Zhang } 174807171c5SHong Zhang c[i] = r1; 175807171c5SHong Zhang c[am + i] = r2; 176807171c5SHong Zhang c[am2 + i] = r3; 177807171c5SHong Zhang c[am3 + i] = r4; 178807171c5SHong Zhang } 179807171c5SHong Zhang b1 += bm4; 180807171c5SHong Zhang b2 += bm4; 181807171c5SHong Zhang b3 += bm4; 182807171c5SHong Zhang b4 += bm4; 183807171c5SHong Zhang 184807171c5SHong Zhang /* RAB[:,col] = R*C[:,col] */ 185807171c5SHong Zhang colrm = col*rm; 186807171c5SHong Zhang for (i=0; i<rm; i++) { /* over rows of R in those columns */ 187807171c5SHong Zhang r1 = r2 = r3 = r4 = 0.0; 188807171c5SHong Zhang n = r->i[i+1] - r->i[i]; 189807171c5SHong Zhang rj = r->j + r->i[i]; 190807171c5SHong Zhang ra = r->a + r->i[i]; 191807171c5SHong Zhang for (j=0; j<n; j++) { 192807171c5SHong Zhang r1 += (*ra)*c[*rj]; 193807171c5SHong Zhang r2 += (*ra)*c2[*rj]; 194807171c5SHong Zhang r3 += (*ra)*c3[*rj]; 195807171c5SHong Zhang r4 += (*ra++)*c4[*rj++]; 196807171c5SHong Zhang } 197807171c5SHong Zhang d[colrm + i] = r1; 198807171c5SHong Zhang d[colrm + rm + i] = r2; 199807171c5SHong Zhang d[colrm + rm2 + i] = r3; 200807171c5SHong Zhang d[colrm + rm3 + i] = r4; 201807171c5SHong Zhang } 202807171c5SHong Zhang } 203807171c5SHong Zhang for (; col<cn; col++) { /* over extra columns of C */ 204807171c5SHong Zhang for (i=0; i<am; i++) { /* over rows of A in those columns */ 205807171c5SHong Zhang r1 = 0.0; 206807171c5SHong Zhang n = a->i[i+1] - a->i[i]; 207807171c5SHong Zhang aj = a->j + a->i[i]; 208807171c5SHong Zhang aa = a->a + a->i[i]; 209807171c5SHong Zhang for (j=0; j<n; j++) { 210807171c5SHong Zhang r1 += (*aa++)*b1[*aj++]; 211807171c5SHong Zhang } 212807171c5SHong Zhang c[i] = r1; 213807171c5SHong Zhang } 214807171c5SHong Zhang b1 += bm; 215807171c5SHong Zhang 216807171c5SHong Zhang for (i=0; i<rm; i++) { /* over rows of R in those columns */ 217807171c5SHong Zhang r1 = 0.0; 218807171c5SHong Zhang n = r->i[i+1] - r->i[i]; 219807171c5SHong Zhang rj = r->j + r->i[i]; 220807171c5SHong Zhang ra = r->a + r->i[i]; 221807171c5SHong Zhang for (j=0; j<n; j++) { 222807171c5SHong Zhang r1 += (*ra++)*c[*rj++]; 223807171c5SHong Zhang } 224807171c5SHong Zhang d[col*rm + i] = r1; 225807171c5SHong Zhang } 226807171c5SHong Zhang } 227807171c5SHong Zhang ierr = PetscLogFlops(cn*2.0*(a->nz + r->nz));CHKERRQ(ierr); 228807171c5SHong Zhang 2298c778c55SBarry Smith ierr = MatDenseRestoreArray(B,&b);CHKERRQ(ierr); 2308c778c55SBarry Smith ierr = MatDenseRestoreArray(RAB,&d);CHKERRQ(ierr); 231807171c5SHong Zhang ierr = MatAssemblyBegin(RAB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 232807171c5SHong Zhang ierr = MatAssemblyEnd(RAB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 233807171c5SHong Zhang PetscFunctionReturn(0); 234807171c5SHong Zhang } 235807171c5SHong Zhang 236807171c5SHong Zhang #undef __FUNCT__ 23755bea0ebSHong Zhang #define __FUNCT__ "MatRARtNumeric_SeqAIJ_SeqAIJ_colorrart" 23855bea0ebSHong Zhang PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ_colorrart(Mat A,Mat R,Mat C) 239807171c5SHong Zhang { 240807171c5SHong Zhang PetscErrorCode ierr; 241257c235dSHong Zhang Mat_SeqAIJ *c = (Mat_SeqAIJ*)C->data; 242257c235dSHong Zhang Mat_RARt *rart=c->rart; 243807171c5SHong Zhang MatTransposeColoring matcoloring; 2448d0a38e4SHong Zhang Mat Rt,RARt; 245807171c5SHong Zhang 246807171c5SHong Zhang PetscFunctionBegin; 247807171c5SHong Zhang /* Get dense Rt by Apply MatTransposeColoring to R */ 248807171c5SHong Zhang matcoloring = rart->matcoloring; 249807171c5SHong Zhang Rt = rart->Rt; 250807171c5SHong Zhang ierr = MatTransColoringApplySpToDen(matcoloring,R,Rt);CHKERRQ(ierr); 251807171c5SHong Zhang 25250647e95SHong Zhang /* Get dense RARt = R*A*Rt -- dominates! */ 253807171c5SHong Zhang RARt = rart->RARt; 254807171c5SHong Zhang ierr = MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense(R,A,Rt,RARt,rart->work);CHKERRQ(ierr); 255807171c5SHong Zhang 256807171c5SHong Zhang /* Recover C from C_dense */ 257807171c5SHong Zhang ierr = MatTransColoringApplyDenToSp(matcoloring,RARt,C);CHKERRQ(ierr); 258807171c5SHong Zhang PetscFunctionReturn(0); 259807171c5SHong Zhang } 260807171c5SHong Zhang 261807171c5SHong Zhang #undef __FUNCT__ 26255bea0ebSHong Zhang #define __FUNCT__ "MatRARtSymbolic_SeqAIJ_SeqAIJ_matmattransposemult" 26355bea0ebSHong Zhang PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ_matmattransposemult(Mat A,Mat R,PetscReal fill,Mat *C) 264807171c5SHong Zhang { 265807171c5SHong Zhang PetscErrorCode ierr; 26655bea0ebSHong Zhang Mat ARt,RARt; 2674fa85fdeSHong Zhang Mat_SeqAIJ *c; 2684fa85fdeSHong Zhang Mat_RARt *rart; 2694fa85fdeSHong Zhang 27095a72cc5SHong Zhang PetscFunctionBegin; 271b3b4fc5aSHong Zhang /* must use '-mat_no_inode' with '-matmattransmult_color 1' - do not knwo why? */ 27236104f73SHong Zhang ierr = MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(A,R,fill,&ARt);CHKERRQ(ierr); 27336104f73SHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(R,ARt,fill,&RARt);CHKERRQ(ierr); 2743b1b9624SHong Zhang *C = RARt; 27555bea0ebSHong Zhang RARt->ops->rartnumeric = MatRARtNumeric_SeqAIJ_SeqAIJ_matmattransposemult; 27655bea0ebSHong Zhang 2773b1b9624SHong Zhang ierr = PetscNew(Mat_RARt,&rart);CHKERRQ(ierr); 2783b1b9624SHong Zhang c = (Mat_SeqAIJ*)(*C)->data; 2793b1b9624SHong Zhang c->rart = rart; 2803b1b9624SHong Zhang rart->ARt = ARt; 2813b1b9624SHong Zhang rart->destroy = RARt->ops->destroy; 2823b1b9624SHong Zhang RARt->ops->destroy = MatDestroy_SeqAIJ_RARt; 2832cfd4408SHong Zhang #if defined(PETSC_USE_INFO) 2842cfd4408SHong 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); 2852cfd4408SHong Zhang #endif 28655bea0ebSHong Zhang PetscFunctionReturn(0); 28736104f73SHong Zhang } 28855bea0ebSHong Zhang 28955bea0ebSHong Zhang #undef __FUNCT__ 29055bea0ebSHong Zhang #define __FUNCT__ "MatRARtNumeric_SeqAIJ_SeqAIJ_matmattransposemult" 29155bea0ebSHong Zhang PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ_matmattransposemult(Mat A,Mat R,Mat C) 29255bea0ebSHong Zhang { 29355bea0ebSHong Zhang PetscErrorCode ierr; 29455bea0ebSHong Zhang Mat_SeqAIJ *c=(Mat_SeqAIJ*)C->data; 29555bea0ebSHong Zhang Mat_RARt *rart=c->rart; 29655bea0ebSHong Zhang Mat ARt=rart->ARt; 29755bea0ebSHong Zhang 29855bea0ebSHong Zhang PetscFunctionBegin; 299b3b4fc5aSHong Zhang ierr = MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(A,R,ARt);CHKERRQ(ierr); /* dominate! */ 30055bea0ebSHong Zhang ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ(R,ARt,C);CHKERRQ(ierr); 30155bea0ebSHong Zhang PetscFunctionReturn(0); 302807171c5SHong Zhang } 30355bea0ebSHong Zhang 30455bea0ebSHong Zhang #undef __FUNCT__ 30555bea0ebSHong Zhang #define __FUNCT__ "MatRARtSymbolic_SeqAIJ_SeqAIJ" 30655bea0ebSHong Zhang PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat R,PetscReal fill,Mat *C) 30755bea0ebSHong Zhang { 30855bea0ebSHong Zhang PetscErrorCode ierr; 30995a72cc5SHong Zhang Mat Rt; 31055bea0ebSHong Zhang Mat_SeqAIJ *c; 31155bea0ebSHong Zhang Mat_RARt *rart; 31255bea0ebSHong Zhang 31355bea0ebSHong Zhang PetscFunctionBegin; 3145008f5a7SHong Zhang ierr = MatTranspose_SeqAIJ(R,MAT_INITIAL_MATRIX,&Rt);CHKERRQ(ierr); 3155008f5a7SHong Zhang ierr = MatMatMatMultSymbolic_SeqAIJ_SeqAIJ_SeqAIJ(R,A,Rt,fill,C);CHKERRQ(ierr); 31695a72cc5SHong Zhang 31795a72cc5SHong Zhang ierr = PetscNew(Mat_RARt,&rart);CHKERRQ(ierr); 31895a72cc5SHong Zhang rart->Rt = Rt; 31995a72cc5SHong Zhang c = (Mat_SeqAIJ*)(*C)->data; 32095a72cc5SHong Zhang c->rart = rart; 32195a72cc5SHong Zhang rart->destroy = (*C)->ops->destroy; 32295a72cc5SHong Zhang (*C)->ops->destroy = MatDestroy_SeqAIJ_RARt; 32355bea0ebSHong Zhang (*C)->ops->rartnumeric = MatRARtNumeric_SeqAIJ_SeqAIJ; 3242cfd4408SHong Zhang #if defined(PETSC_USE_INFO) 3252cfd4408SHong Zhang ierr = PetscInfo(*C,"Use Rt=R^T and C=R*A*Rt via MatMatMatMult() to avoid sparse inner products\n");CHKERRQ(ierr); 3262cfd4408SHong Zhang #endif 32755bea0ebSHong Zhang PetscFunctionReturn(0); 32855bea0ebSHong Zhang } 32955bea0ebSHong Zhang 33055bea0ebSHong Zhang #undef __FUNCT__ 33155bea0ebSHong Zhang #define __FUNCT__ "MatRARtNumeric_SeqAIJ_SeqAIJ" 33255bea0ebSHong Zhang PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ(Mat A,Mat R,Mat C) 33355bea0ebSHong Zhang { 33455bea0ebSHong Zhang PetscErrorCode ierr; 33555bea0ebSHong Zhang Mat_SeqAIJ *c = (Mat_SeqAIJ*)C->data; 33655bea0ebSHong Zhang Mat_RARt *rart = c->rart; 33755bea0ebSHong Zhang Mat Rt = rart->Rt; 33855bea0ebSHong Zhang 33955bea0ebSHong Zhang PetscFunctionBegin; 34055bea0ebSHong Zhang ierr = MatTranspose_SeqAIJ(R,MAT_REUSE_MATRIX,&Rt);CHKERRQ(ierr); 34155bea0ebSHong Zhang ierr = MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqAIJ(R,A,Rt,C);CHKERRQ(ierr); 34255bea0ebSHong Zhang PetscFunctionReturn(0); 34355bea0ebSHong Zhang } 34455bea0ebSHong Zhang 34555bea0ebSHong Zhang #undef __FUNCT__ 34655bea0ebSHong Zhang #define __FUNCT__ "MatRARt_SeqAIJ_SeqAIJ" 34755bea0ebSHong Zhang PetscErrorCode MatRARt_SeqAIJ_SeqAIJ(Mat A,Mat R,MatReuse scall,PetscReal fill,Mat *C) 34855bea0ebSHong Zhang { 34955bea0ebSHong Zhang PetscErrorCode ierr; 35055bea0ebSHong Zhang const char *algTypes[3] = {"matmatmatmult","matmattransposemult","coloring_rart"}; 35155bea0ebSHong Zhang PetscInt alg=0; /* set default algorithm */ 35255bea0ebSHong Zhang 35355bea0ebSHong Zhang PetscFunctionBegin; 35455bea0ebSHong Zhang if (scall == MAT_INITIAL_MATRIX) { 355b56132cbSHong Zhang ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr); 35655bea0ebSHong Zhang ierr = PetscOptionsEList("-matrart_via","Algorithmic approach","MatRARt",algTypes,3,algTypes[0],&alg,NULL);CHKERRQ(ierr); 357b56132cbSHong Zhang ierr = PetscOptionsEnd();CHKERRQ(ierr); 358b56132cbSHong Zhang 35955bea0ebSHong Zhang ierr = PetscLogEventBegin(MAT_RARtSymbolic,A,R,0,0);CHKERRQ(ierr); 36055bea0ebSHong Zhang switch (alg) { 36155bea0ebSHong Zhang case 1: 36255bea0ebSHong Zhang /* via matmattransposemult: ARt=A*R^T, C=R*ARt - matrix coloring can be applied to A*R^T */ 36355bea0ebSHong Zhang ierr = MatRARtSymbolic_SeqAIJ_SeqAIJ_matmattransposemult(A,R,fill,C);CHKERRQ(ierr); 36455bea0ebSHong Zhang break; 36555bea0ebSHong Zhang case 2: 36655bea0ebSHong Zhang /* via coloring_rart: apply coloring C = R*A*R^T */ 36755bea0ebSHong Zhang ierr = MatRARtSymbolic_SeqAIJ_SeqAIJ_colorrart(A,R,fill,C);CHKERRQ(ierr); 36855bea0ebSHong Zhang break; 36955bea0ebSHong Zhang default: 37055bea0ebSHong Zhang /* via matmatmatmult: Rt=R^T, C=R*A*Rt - avoid inefficient sparse inner products */ 37155bea0ebSHong Zhang ierr = MatRARtSymbolic_SeqAIJ_SeqAIJ(A,R,fill,C);CHKERRQ(ierr); 37255bea0ebSHong Zhang break; 37355bea0ebSHong Zhang } 3745008f5a7SHong Zhang ierr = PetscLogEventEnd(MAT_RARtSymbolic,A,R,0,0);CHKERRQ(ierr); 3755008f5a7SHong Zhang } 3765008f5a7SHong Zhang 3775008f5a7SHong Zhang ierr = PetscLogEventBegin(MAT_RARtNumeric,A,R,0,0);CHKERRQ(ierr); 37855bea0ebSHong Zhang ierr = (*(*C)->ops->rartnumeric)(A,R,*C);CHKERRQ(ierr); 3795008f5a7SHong Zhang ierr = PetscLogEventEnd(MAT_RARtNumeric,A,R,0,0);CHKERRQ(ierr); 380807171c5SHong Zhang PetscFunctionReturn(0); 381807171c5SHong Zhang } 382