1807171c5SHong Zhang 2807171c5SHong Zhang /* 3807171c5SHong Zhang Defines matrix-matrix product routines for pairs of SeqAIJ matrices 4807171c5SHong Zhang C = P * A * P^T 5807171c5SHong Zhang */ 6807171c5SHong 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 /* 13807171c5SHong Zhang MatApplyPAPt_Symbolic_SeqAIJ_SeqAIJ - Forms the symbolic product of two SeqAIJ matrices 14807171c5SHong Zhang C = P * A * P^T; 15807171c5SHong Zhang 16807171c5SHong Zhang Note: C is assumed to be uncreated. 17807171c5SHong Zhang If this is not the case, Destroy C before calling this routine. 18807171c5SHong Zhang */ 19807171c5SHong Zhang #undef __FUNCT__ 20807171c5SHong Zhang #define __FUNCT__ "MatApplyPAPt_Symbolic_SeqAIJ_SeqAIJ" 21807171c5SHong Zhang PetscErrorCode MatApplyPAPt_Symbolic_SeqAIJ_SeqAIJ(Mat A,Mat P,Mat *C) 22807171c5SHong Zhang { 23807171c5SHong Zhang /* Note: This code is virtually identical to that of MatApplyPtAP_SeqAIJ_Symbolic */ 24807171c5SHong Zhang /* and MatMatMult_SeqAIJ_SeqAIJ_Symbolic. Perhaps they could be merged nicely. */ 25807171c5SHong Zhang PetscErrorCode ierr; 260298fd71SBarry Smith PetscFreeSpaceList free_space=NULL,current_space=NULL; 27807171c5SHong Zhang Mat_SeqAIJ *a =(Mat_SeqAIJ*)A->data,*p=(Mat_SeqAIJ*)P->data,*c; 28807171c5SHong Zhang PetscInt *ai =a->i,*aj=a->j,*ajj,*pi=p->i,*pj=p->j,*pti,*ptj,*ptjj; 29807171c5SHong Zhang PetscInt *ci,*cj,*paj,*padenserow,*pasparserow,*denserow,*sparserow; 30807171c5SHong Zhang PetscInt an=A->cmap->N,am=A->rmap->N,pn=P->cmap->N,pm=P->rmap->N; 31807171c5SHong Zhang PetscInt i,j,k,pnzi,arow,anzj,panzi,ptrow,ptnzj,cnzi; 32807171c5SHong Zhang MatScalar *ca; 33807171c5SHong Zhang 34807171c5SHong Zhang PetscFunctionBegin; 35807171c5SHong Zhang /* some error checking which could be moved into interface layer */ 36807171c5SHong Zhang if (pn!=am) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %D != %D",pn,am); 37807171c5SHong Zhang if (am!=an) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Matrix 'A' must be square, %D != %D",am, an); 38807171c5SHong Zhang 39807171c5SHong Zhang /* Set up timers */ 40807171c5SHong Zhang ierr = PetscLogEventBegin(MAT_Applypapt_symbolic,A,P,0,0);CHKERRQ(ierr); 41807171c5SHong Zhang 42807171c5SHong Zhang /* Create ij structure of P^T */ 43807171c5SHong Zhang ierr = MatGetSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr); 44807171c5SHong Zhang 45807171c5SHong Zhang /* Allocate ci array, arrays for fill computation and */ 46807171c5SHong Zhang /* free space for accumulating nonzero column info */ 47807171c5SHong Zhang ierr = PetscMalloc(((pm+1)*1)*sizeof(PetscInt),&ci);CHKERRQ(ierr); 48807171c5SHong Zhang ci[0] = 0; 49807171c5SHong Zhang 50807171c5SHong Zhang ierr = PetscMalloc4(an,PetscInt,&padenserow,an,PetscInt,&pasparserow,pm,PetscInt,&denserow,pm,PetscInt,&sparserow);CHKERRQ(ierr); 51807171c5SHong Zhang ierr = PetscMemzero(padenserow,an*sizeof(PetscInt));CHKERRQ(ierr); 52807171c5SHong Zhang ierr = PetscMemzero(pasparserow,an*sizeof(PetscInt));CHKERRQ(ierr); 53807171c5SHong Zhang ierr = PetscMemzero(denserow,pm*sizeof(PetscInt));CHKERRQ(ierr); 54807171c5SHong Zhang ierr = PetscMemzero(sparserow,pm*sizeof(PetscInt));CHKERRQ(ierr); 55807171c5SHong Zhang 56807171c5SHong Zhang /* Set initial free space to be nnz(A) scaled by aspect ratio of Pt. */ 57807171c5SHong Zhang /* This should be reasonable if sparsity of PAPt is similar to that of A. */ 58a2ea699eSBarry Smith ierr = PetscFreeSpaceGet((ai[am]/pn)*pm,&free_space);CHKERRQ(ierr); 59807171c5SHong Zhang current_space = free_space; 60807171c5SHong Zhang 61807171c5SHong Zhang /* Determine fill for each row of C: */ 62807171c5SHong Zhang for (i=0; i<pm; i++) { 63807171c5SHong Zhang pnzi = pi[i+1] - pi[i]; 64807171c5SHong Zhang panzi = 0; 65807171c5SHong Zhang /* Get symbolic sparse row of PA: */ 66807171c5SHong Zhang for (j=0; j<pnzi; j++) { 67807171c5SHong Zhang arow = *pj++; 68807171c5SHong Zhang anzj = ai[arow+1] - ai[arow]; 69807171c5SHong Zhang ajj = aj + ai[arow]; 70807171c5SHong Zhang for (k=0; k<anzj; k++) { 71807171c5SHong Zhang if (!padenserow[ajj[k]]) { 72807171c5SHong Zhang padenserow[ajj[k]] = -1; 73807171c5SHong Zhang pasparserow[panzi++] = ajj[k]; 74807171c5SHong Zhang } 75807171c5SHong Zhang } 76807171c5SHong Zhang } 77807171c5SHong Zhang /* Using symbolic row of PA, determine symbolic row of C: */ 78807171c5SHong Zhang paj = pasparserow; 79807171c5SHong Zhang cnzi = 0; 80807171c5SHong Zhang for (j=0; j<panzi; j++) { 81807171c5SHong Zhang ptrow = *paj++; 82807171c5SHong Zhang ptnzj = pti[ptrow+1] - pti[ptrow]; 83807171c5SHong Zhang ptjj = ptj + pti[ptrow]; 84807171c5SHong Zhang for (k=0; k<ptnzj; k++) { 85807171c5SHong Zhang if (!denserow[ptjj[k]]) { 86807171c5SHong Zhang denserow[ptjj[k]] = -1; 87807171c5SHong Zhang sparserow[cnzi++] = ptjj[k]; 88807171c5SHong Zhang } 89807171c5SHong Zhang } 90807171c5SHong Zhang } 91807171c5SHong Zhang 92807171c5SHong Zhang /* sort sparse representation */ 93807171c5SHong Zhang ierr = PetscSortInt(cnzi,sparserow);CHKERRQ(ierr); 94807171c5SHong Zhang 95807171c5SHong Zhang /* If free space is not available, make more free space */ 96807171c5SHong Zhang /* Double the amount of total space in the list */ 97807171c5SHong Zhang if (current_space->local_remaining<cnzi) { 98807171c5SHong Zhang ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 99807171c5SHong Zhang } 100807171c5SHong Zhang 101807171c5SHong Zhang /* Copy data into free space, and zero out dense row */ 102807171c5SHong Zhang ierr = PetscMemcpy(current_space->array,sparserow,cnzi*sizeof(PetscInt));CHKERRQ(ierr); 1032205254eSKarl Rupp 104807171c5SHong Zhang current_space->array += cnzi; 105807171c5SHong Zhang current_space->local_used += cnzi; 106807171c5SHong Zhang current_space->local_remaining -= cnzi; 107807171c5SHong Zhang 108807171c5SHong Zhang for (j=0; j<panzi; j++) { 109807171c5SHong Zhang padenserow[pasparserow[j]] = 0; 110807171c5SHong Zhang } 111807171c5SHong Zhang for (j=0; j<cnzi; j++) { 112807171c5SHong Zhang denserow[sparserow[j]] = 0; 113807171c5SHong Zhang } 114807171c5SHong Zhang ci[i+1] = ci[i] + cnzi; 115807171c5SHong Zhang } 116807171c5SHong Zhang /* column indices are in the list of free space */ 117807171c5SHong Zhang /* Allocate space for cj, initialize cj, and */ 118807171c5SHong Zhang /* destroy list of free space and other temporary array(s) */ 119807171c5SHong Zhang ierr = PetscMalloc((ci[pm]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr); 120807171c5SHong Zhang ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr); 121807171c5SHong Zhang ierr = PetscFree4(padenserow,pasparserow,denserow,sparserow);CHKERRQ(ierr); 122807171c5SHong Zhang 123807171c5SHong Zhang /* Allocate space for ca */ 124807171c5SHong Zhang ierr = PetscMalloc((ci[pm]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr); 125807171c5SHong Zhang ierr = PetscMemzero(ca,(ci[pm]+1)*sizeof(MatScalar));CHKERRQ(ierr); 126807171c5SHong Zhang 127807171c5SHong Zhang /* put together the new matrix */ 128ce94432eSBarry Smith ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),pm,pm,ci,cj,ca,C);CHKERRQ(ierr); 129a2f3521dSMark F. Adams (*C)->rmap->bs = P->cmap->bs; 130a2f3521dSMark F. Adams (*C)->cmap->bs = P->cmap->bs; 131807171c5SHong Zhang 132807171c5SHong Zhang /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 133807171c5SHong Zhang /* Since these are PETSc arrays, change flags to free them as necessary. */ 134807171c5SHong Zhang c = (Mat_SeqAIJ*)((*C)->data); 135807171c5SHong Zhang c->free_a = PETSC_TRUE; 136807171c5SHong Zhang c->free_ij = PETSC_TRUE; 137807171c5SHong Zhang c->nonew = 0; 138807171c5SHong Zhang 139807171c5SHong Zhang /* Clean up. */ 140807171c5SHong Zhang ierr = MatRestoreSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr); 141807171c5SHong Zhang 142807171c5SHong Zhang ierr = PetscLogEventEnd(MAT_Applypapt_symbolic,A,P,0,0);CHKERRQ(ierr); 143807171c5SHong Zhang PetscFunctionReturn(0); 144807171c5SHong Zhang } 145807171c5SHong Zhang 146807171c5SHong Zhang /* 147807171c5SHong Zhang MatApplyPAPt_Numeric_SeqAIJ - Forms the numeric product of two SeqAIJ matrices 148807171c5SHong Zhang C = P * A * P^T; 149807171c5SHong Zhang Note: C must have been created by calling MatApplyPAPt_Symbolic_SeqAIJ. 150807171c5SHong Zhang */ 151807171c5SHong Zhang #undef __FUNCT__ 152807171c5SHong Zhang #define __FUNCT__ "MatApplyPAPt_Numeric_SeqAIJ_SeqAIJ" 153807171c5SHong Zhang PetscErrorCode MatApplyPAPt_Numeric_SeqAIJ_SeqAIJ(Mat A,Mat P,Mat C) 154807171c5SHong Zhang { 155807171c5SHong Zhang PetscErrorCode ierr; 156807171c5SHong Zhang PetscInt flops=0; 157807171c5SHong Zhang Mat_SeqAIJ *a = (Mat_SeqAIJ*) A->data; 158807171c5SHong Zhang Mat_SeqAIJ *p = (Mat_SeqAIJ*) P->data; 159807171c5SHong Zhang Mat_SeqAIJ *c = (Mat_SeqAIJ*) C->data; 160807171c5SHong Zhang PetscInt *ai = a->i,*aj=a->j,*ajj,*pi=p->i,*pj=p->j,*pjj=p->j,*paj,*pajdense,*ptj; 161807171c5SHong Zhang PetscInt *ci = c->i,*cj=c->j; 162807171c5SHong Zhang PetscInt an = A->cmap->N,am=A->rmap->N,pn=P->cmap->N,pm=P->rmap->N,cn=C->cmap->N,cm=C->rmap->N; 163807171c5SHong Zhang PetscInt i,j,k,k1,k2,pnzi,anzj,panzj,arow,ptcol,ptnzj,cnzi; 164807171c5SHong Zhang MatScalar *aa=a->a,*pa=p->a,*pta=p->a,*ptaj,*paa,*aaj,*ca=c->a,sum; 165807171c5SHong Zhang 166807171c5SHong Zhang PetscFunctionBegin; 167807171c5SHong Zhang /* This error checking should be unnecessary if the symbolic was performed */ 168807171c5SHong Zhang if (pm!=cm) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %D != %D",pm,cm); 169807171c5SHong Zhang if (pn!=am) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %D != %D",pn,am); 170807171c5SHong Zhang if (am!=an) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Matrix 'A' must be square, %D != %D",am, an); 171807171c5SHong Zhang if (pm!=cn) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %D != %D",pm, cn); 172807171c5SHong Zhang 173807171c5SHong Zhang /* Set up timers */ 174807171c5SHong Zhang ierr = PetscLogEventBegin(MAT_Applypapt_numeric,A,P,C,0);CHKERRQ(ierr); 175807171c5SHong Zhang ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr); 176807171c5SHong Zhang 177807171c5SHong Zhang ierr = PetscMalloc3(an,MatScalar,&paa,an,PetscInt,&paj,an,PetscInt,&pajdense);CHKERRQ(ierr); 178807171c5SHong Zhang ierr = PetscMemzero(paa,an*(sizeof(MatScalar)+2*sizeof(PetscInt)));CHKERRQ(ierr); 179807171c5SHong Zhang 180807171c5SHong Zhang for (i=0; i<pm; i++) { 181807171c5SHong Zhang /* Form sparse row of P*A */ 182807171c5SHong Zhang pnzi = pi[i+1] - pi[i]; 183807171c5SHong Zhang panzj = 0; 184807171c5SHong Zhang for (j=0; j<pnzi; j++) { 185807171c5SHong Zhang arow = *pj++; 186807171c5SHong Zhang anzj = ai[arow+1] - ai[arow]; 187807171c5SHong Zhang ajj = aj + ai[arow]; 188807171c5SHong Zhang aaj = aa + ai[arow]; 189807171c5SHong Zhang for (k=0; k<anzj; k++) { 190807171c5SHong Zhang if (!pajdense[ajj[k]]) { 191807171c5SHong Zhang pajdense[ajj[k]] = -1; 192807171c5SHong Zhang paj[panzj++] = ajj[k]; 193807171c5SHong Zhang } 194807171c5SHong Zhang paa[ajj[k]] += (*pa)*aaj[k]; 195807171c5SHong Zhang } 196807171c5SHong Zhang flops += 2*anzj; 197807171c5SHong Zhang pa++; 198807171c5SHong Zhang } 199807171c5SHong Zhang 200807171c5SHong Zhang /* Sort the j index array for quick sparse axpy. */ 201807171c5SHong Zhang ierr = PetscSortInt(panzj,paj);CHKERRQ(ierr); 202807171c5SHong Zhang 203807171c5SHong Zhang /* Compute P*A*P^T using sparse inner products. */ 204807171c5SHong Zhang /* Take advantage of pre-computed (i,j) of C for locations of non-zeros. */ 205807171c5SHong Zhang cnzi = ci[i+1] - ci[i]; 206807171c5SHong Zhang for (j=0; j<cnzi; j++) { 207807171c5SHong Zhang /* Form sparse inner product of current row of P*A with (*cj++) col of P^T. */ 208807171c5SHong Zhang ptcol = *cj++; 209807171c5SHong Zhang ptnzj = pi[ptcol+1] - pi[ptcol]; 210807171c5SHong Zhang ptj = pjj + pi[ptcol]; 211807171c5SHong Zhang ptaj = pta + pi[ptcol]; 212807171c5SHong Zhang sum = 0.; 213807171c5SHong Zhang k1 = 0; 214807171c5SHong Zhang k2 = 0; 215807171c5SHong Zhang while ((k1<panzj) && (k2<ptnzj)) { 216807171c5SHong Zhang if (paj[k1]==ptj[k2]) { 217807171c5SHong Zhang sum += paa[paj[k1++]]*ptaj[k2++]; 218807171c5SHong Zhang } else if (paj[k1] < ptj[k2]) { 219807171c5SHong Zhang k1++; 220807171c5SHong Zhang } else /* if (paj[k1] > ptj[k2]) */ { 221807171c5SHong Zhang k2++; 222807171c5SHong Zhang } 223807171c5SHong Zhang } 224807171c5SHong Zhang *ca++ = sum; 225807171c5SHong Zhang } 226807171c5SHong Zhang 227807171c5SHong Zhang /* Zero the current row info for P*A */ 228807171c5SHong Zhang for (j=0;j<panzj;j++) { 229807171c5SHong Zhang paa[paj[j]] = 0.; 230807171c5SHong Zhang pajdense[paj[j]] = 0; 231807171c5SHong Zhang } 232807171c5SHong Zhang } 233807171c5SHong Zhang 234807171c5SHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 235807171c5SHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 236807171c5SHong Zhang ierr = PetscFree3(paa,paj,pajdense);CHKERRQ(ierr); 237807171c5SHong Zhang ierr = PetscLogFlops(flops);CHKERRQ(ierr); 238807171c5SHong Zhang ierr = PetscLogEventEnd(MAT_Applypapt_numeric,A,P,C,0);CHKERRQ(ierr); 239807171c5SHong Zhang PetscFunctionReturn(0); 240807171c5SHong Zhang } 241807171c5SHong Zhang 242807171c5SHong Zhang #undef __FUNCT__ 243807171c5SHong Zhang #define __FUNCT__ "MatApplyPAPt_SeqAIJ_SeqAIJ" 244807171c5SHong Zhang PetscErrorCode MatApplyPAPt_SeqAIJ_SeqAIJ(Mat A,Mat P,Mat *C) 245807171c5SHong Zhang { 246807171c5SHong Zhang PetscErrorCode ierr; 247807171c5SHong Zhang 248807171c5SHong Zhang PetscFunctionBegin; 249807171c5SHong Zhang ierr = PetscLogEventBegin(MAT_Applypapt,A,P,0,0);CHKERRQ(ierr); 250807171c5SHong Zhang ierr = MatApplyPAPt_Symbolic_SeqAIJ_SeqAIJ(A,P,C);CHKERRQ(ierr); 251807171c5SHong Zhang ierr = MatApplyPAPt_Numeric_SeqAIJ_SeqAIJ(A,P,*C);CHKERRQ(ierr); 252807171c5SHong Zhang ierr = PetscLogEventEnd(MAT_Applypapt,A,P,0,0);CHKERRQ(ierr); 253807171c5SHong Zhang PetscFunctionReturn(0); 254807171c5SHong Zhang } 255807171c5SHong Zhang 256807171c5SHong Zhang /*--------------------------------------------------*/ 257807171c5SHong Zhang /* 258807171c5SHong Zhang Defines projective product routines where A is a SeqAIJ matrix 259807171c5SHong Zhang C = R * A * R^T 260807171c5SHong Zhang */ 261807171c5SHong Zhang 262807171c5SHong Zhang #undef __FUNCT__ 263807171c5SHong Zhang #define __FUNCT__ "PetscContainerDestroy_Mat_RARt" 264807171c5SHong Zhang PetscErrorCode PetscContainerDestroy_Mat_RARt(void *ptr) 265807171c5SHong Zhang { 266807171c5SHong Zhang PetscErrorCode ierr; 267807171c5SHong Zhang Mat_RARt *rart=(Mat_RARt*)ptr; 268807171c5SHong Zhang 269807171c5SHong Zhang PetscFunctionBegin; 270807171c5SHong Zhang ierr = MatTransposeColoringDestroy(&rart->matcoloring);CHKERRQ(ierr); 271807171c5SHong Zhang ierr = MatDestroy(&rart->Rt);CHKERRQ(ierr); 272807171c5SHong Zhang ierr = MatDestroy(&rart->RARt);CHKERRQ(ierr); 273807171c5SHong Zhang ierr = PetscFree(rart->work);CHKERRQ(ierr); 274807171c5SHong Zhang ierr = PetscFree(rart);CHKERRQ(ierr); 275807171c5SHong Zhang PetscFunctionReturn(0); 276807171c5SHong Zhang } 277807171c5SHong Zhang 278807171c5SHong Zhang #undef __FUNCT__ 279807171c5SHong Zhang #define __FUNCT__ "MatDestroy_SeqAIJ_RARt" 280807171c5SHong Zhang PetscErrorCode MatDestroy_SeqAIJ_RARt(Mat A) 281807171c5SHong Zhang { 282807171c5SHong Zhang PetscErrorCode ierr; 283807171c5SHong Zhang PetscContainer container; 2840298fd71SBarry Smith Mat_RARt *rart=NULL; 285807171c5SHong Zhang 286807171c5SHong Zhang PetscFunctionBegin; 287807171c5SHong Zhang ierr = PetscObjectQuery((PetscObject)A,"Mat_RARt",(PetscObject*)&container);CHKERRQ(ierr); 288807171c5SHong Zhang if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit"); 289807171c5SHong Zhang ierr = PetscContainerGetPointer(container,(void**)&rart);CHKERRQ(ierr); 290807171c5SHong Zhang A->ops->destroy = rart->destroy; 291807171c5SHong Zhang if (A->ops->destroy) { 292807171c5SHong Zhang ierr = (*A->ops->destroy)(A);CHKERRQ(ierr); 293807171c5SHong Zhang } 294807171c5SHong Zhang ierr = PetscObjectCompose((PetscObject)A,"Mat_RARt",0);CHKERRQ(ierr); 295807171c5SHong Zhang PetscFunctionReturn(0); 296807171c5SHong Zhang } 297807171c5SHong Zhang 298807171c5SHong Zhang #undef __FUNCT__ 299807171c5SHong Zhang #define __FUNCT__ "MatRARtSymbolic_SeqAIJ_SeqAIJ" 300807171c5SHong Zhang PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat R,PetscReal fill,Mat *C) 301807171c5SHong Zhang { 302807171c5SHong Zhang PetscErrorCode ierr; 303807171c5SHong Zhang Mat P; 304807171c5SHong Zhang PetscInt *rti,*rtj; 305807171c5SHong Zhang Mat_RARt *rart; 306807171c5SHong Zhang PetscContainer container; 307807171c5SHong Zhang MatTransposeColoring matcoloring; 308807171c5SHong Zhang ISColoring iscoloring; 3098d0a38e4SHong Zhang Mat Rt_dense,RARt_dense; 310807171c5SHong Zhang PetscLogDouble GColor=0.0,MCCreate=0.0,MDenCreate=0.0,t0,tf,etime=0.0; 311807171c5SHong Zhang Mat_SeqAIJ *c; 312807171c5SHong Zhang 313807171c5SHong Zhang PetscFunctionBegin; 3148563dfccSBarry Smith ierr = PetscTime(&t0);CHKERRQ(ierr); 315807171c5SHong Zhang /* create symbolic P=Rt */ 316807171c5SHong Zhang ierr = MatGetSymbolicTranspose_SeqAIJ(R,&rti,&rtj);CHKERRQ(ierr); 3170298fd71SBarry Smith ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,R->cmap->n,R->rmap->n,rti,rtj,NULL,&P);CHKERRQ(ierr); 318807171c5SHong Zhang 319807171c5SHong Zhang /* get symbolic C=Pt*A*P */ 320807171c5SHong Zhang ierr = MatPtAPSymbolic_SeqAIJ_SeqAIJ(A,P,fill,C);CHKERRQ(ierr); 321fd8e3dafSHong Zhang (*C)->rmap->bs = R->rmap->bs; 322fd8e3dafSHong Zhang (*C)->cmap->bs = R->rmap->bs; 323807171c5SHong Zhang 324807171c5SHong Zhang /* create a supporting struct */ 325807171c5SHong Zhang ierr = PetscNew(Mat_RARt,&rart);CHKERRQ(ierr); 326807171c5SHong Zhang 327807171c5SHong Zhang /* attach the supporting struct to C */ 328807171c5SHong Zhang ierr = PetscContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr); 329807171c5SHong Zhang ierr = PetscContainerSetPointer(container,rart);CHKERRQ(ierr); 330807171c5SHong Zhang ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_Mat_RARt);CHKERRQ(ierr); 331807171c5SHong Zhang ierr = PetscObjectCompose((PetscObject)(*C),"Mat_RARt",(PetscObject)container);CHKERRQ(ierr); 332807171c5SHong Zhang ierr = PetscContainerDestroy(&container);CHKERRQ(ierr); 3338563dfccSBarry Smith ierr = PetscTime(&tf);CHKERRQ(ierr); 334807171c5SHong Zhang etime += tf - t0; 335807171c5SHong Zhang 336*22e94b5dSHong Zhang /* ------ Use coloring ---------- */ 337807171c5SHong Zhang /* Create MatTransposeColoring from symbolic C=R*A*R^T */ 338807171c5SHong Zhang c = (Mat_SeqAIJ*)(*C)->data; 3398563dfccSBarry Smith ierr = PetscTime(&t0);CHKERRQ(ierr); 340807171c5SHong Zhang ierr = MatGetColoring(*C,MATCOLORINGLF,&iscoloring);CHKERRQ(ierr); 3418563dfccSBarry Smith ierr = PetscTime(&tf);CHKERRQ(ierr); 342807171c5SHong Zhang GColor += tf - t0; 343807171c5SHong Zhang 3448563dfccSBarry Smith ierr = PetscTime(&t0);CHKERRQ(ierr); 345807171c5SHong Zhang ierr = MatTransposeColoringCreate(*C,iscoloring,&matcoloring);CHKERRQ(ierr); 3462205254eSKarl Rupp 347807171c5SHong Zhang rart->matcoloring = matcoloring; 3482205254eSKarl Rupp 349807171c5SHong Zhang ierr = ISColoringDestroy(&iscoloring);CHKERRQ(ierr); 3508563dfccSBarry Smith ierr = PetscTime(&tf);CHKERRQ(ierr); 351807171c5SHong Zhang MCCreate += tf - t0; 352807171c5SHong Zhang 3538563dfccSBarry Smith ierr = PetscTime(&t0);CHKERRQ(ierr); 354807171c5SHong Zhang /* Create Rt_dense */ 355807171c5SHong Zhang ierr = MatCreate(PETSC_COMM_SELF,&Rt_dense);CHKERRQ(ierr); 356807171c5SHong Zhang ierr = MatSetSizes(Rt_dense,A->cmap->n,matcoloring->ncolors,A->cmap->n,matcoloring->ncolors);CHKERRQ(ierr); 357807171c5SHong Zhang ierr = MatSetType(Rt_dense,MATSEQDENSE);CHKERRQ(ierr); 3580298fd71SBarry Smith ierr = MatSeqDenseSetPreallocation(Rt_dense,NULL);CHKERRQ(ierr); 3592205254eSKarl Rupp 360807171c5SHong Zhang Rt_dense->assembled = PETSC_TRUE; 361807171c5SHong Zhang rart->Rt = Rt_dense; 362807171c5SHong Zhang 363807171c5SHong Zhang /* Create RARt_dense = R*A*Rt_dense */ 364807171c5SHong Zhang ierr = MatCreate(PETSC_COMM_SELF,&RARt_dense);CHKERRQ(ierr); 365807171c5SHong Zhang ierr = MatSetSizes(RARt_dense,(*C)->rmap->n,matcoloring->ncolors,(*C)->rmap->n,matcoloring->ncolors);CHKERRQ(ierr); 366807171c5SHong Zhang ierr = MatSetType(RARt_dense,MATSEQDENSE);CHKERRQ(ierr); 3670298fd71SBarry Smith ierr = MatSeqDenseSetPreallocation(RARt_dense,NULL);CHKERRQ(ierr); 3682205254eSKarl Rupp 369807171c5SHong Zhang rart->RARt = RARt_dense; 370807171c5SHong Zhang 371807171c5SHong Zhang /* Allocate work array to store columns of A*R^T used in MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense() */ 372807171c5SHong Zhang ierr = PetscMalloc(A->rmap->n*4*sizeof(PetscScalar),&rart->work);CHKERRQ(ierr); 373807171c5SHong Zhang 3748563dfccSBarry Smith ierr = PetscTime(&tf);CHKERRQ(ierr); 375807171c5SHong Zhang MDenCreate += tf - t0; 376807171c5SHong Zhang 377807171c5SHong Zhang rart->destroy = (*C)->ops->destroy; 378807171c5SHong Zhang (*C)->ops->destroy = MatDestroy_SeqAIJ_RARt; 379807171c5SHong Zhang 380807171c5SHong Zhang /* clean up */ 381807171c5SHong Zhang ierr = MatRestoreSymbolicTranspose_SeqAIJ(R,&rti,&rtj);CHKERRQ(ierr); 382807171c5SHong Zhang ierr = MatDestroy(&P);CHKERRQ(ierr); 383807171c5SHong Zhang 384807171c5SHong Zhang #if defined(PETSC_USE_INFO) 3851ce71dffSSatish Balay { 386807171c5SHong Zhang PetscReal density= (PetscReal)(c->nz)/(RARt_dense->rmap->n*RARt_dense->cmap->n); 387a2ea699eSBarry 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); 388a2ea699eSBarry 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); 3891ce71dffSSatish Balay } 390807171c5SHong Zhang #endif 391807171c5SHong Zhang PetscFunctionReturn(0); 392807171c5SHong Zhang } 393807171c5SHong Zhang 394807171c5SHong Zhang /* 395807171c5SHong Zhang RAB = R * A * B, R and A in seqaij format, B in dense format; 396807171c5SHong Zhang */ 397807171c5SHong Zhang #undef __FUNCT__ 398807171c5SHong Zhang #define __FUNCT__ "MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense" 399807171c5SHong Zhang PetscErrorCode MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense(Mat R,Mat A,Mat B,Mat RAB,PetscScalar *work) 400807171c5SHong Zhang { 401807171c5SHong Zhang Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*r=(Mat_SeqAIJ*)R->data; 402807171c5SHong Zhang PetscErrorCode ierr; 403807171c5SHong Zhang PetscScalar *b,r1,r2,r3,r4,*b1,*b2,*b3,*b4; 404807171c5SHong Zhang MatScalar *aa,*ra; 405807171c5SHong Zhang PetscInt cn =B->cmap->n,bm=B->rmap->n,col,i,j,n,*ai=a->i,*aj,am=A->rmap->n; 406807171c5SHong Zhang PetscInt am2=2*am,am3=3*am,bm4=4*bm; 407807171c5SHong Zhang PetscScalar *d,*c,*c2,*c3,*c4; 408807171c5SHong Zhang PetscInt *rj,rm=R->rmap->n,dm=RAB->rmap->n,dn=RAB->cmap->n; 409807171c5SHong Zhang PetscInt rm2=2*rm,rm3=3*rm,colrm; 410807171c5SHong Zhang 411807171c5SHong Zhang PetscFunctionBegin; 412807171c5SHong Zhang if (!dm || !dn) PetscFunctionReturn(0); 413807171c5SHong 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); 414807171c5SHong 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); 415807171c5SHong 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); 416807171c5SHong 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); 417807171c5SHong Zhang 4188c778c55SBarry Smith ierr = MatDenseGetArray(B,&b);CHKERRQ(ierr); 4198c778c55SBarry Smith ierr = MatDenseGetArray(RAB,&d);CHKERRQ(ierr); 420807171c5SHong Zhang b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm; 421807171c5SHong Zhang c = work; c2 = c + am; c3 = c2 + am; c4 = c3 + am; 422807171c5SHong Zhang for (col=0; col<cn-4; col += 4) { /* over columns of C */ 423807171c5SHong Zhang for (i=0; i<am; i++) { /* over rows of A in those columns */ 424807171c5SHong Zhang r1 = r2 = r3 = r4 = 0.0; 425807171c5SHong Zhang n = ai[i+1] - ai[i]; 426807171c5SHong Zhang aj = a->j + ai[i]; 427807171c5SHong Zhang aa = a->a + ai[i]; 428807171c5SHong Zhang for (j=0; j<n; j++) { 429807171c5SHong Zhang r1 += (*aa)*b1[*aj]; 430807171c5SHong Zhang r2 += (*aa)*b2[*aj]; 431807171c5SHong Zhang r3 += (*aa)*b3[*aj]; 432807171c5SHong Zhang r4 += (*aa++)*b4[*aj++]; 433807171c5SHong Zhang } 434807171c5SHong Zhang c[i] = r1; 435807171c5SHong Zhang c[am + i] = r2; 436807171c5SHong Zhang c[am2 + i] = r3; 437807171c5SHong Zhang c[am3 + i] = r4; 438807171c5SHong Zhang } 439807171c5SHong Zhang b1 += bm4; 440807171c5SHong Zhang b2 += bm4; 441807171c5SHong Zhang b3 += bm4; 442807171c5SHong Zhang b4 += bm4; 443807171c5SHong Zhang 444807171c5SHong Zhang /* RAB[:,col] = R*C[:,col] */ 445807171c5SHong Zhang colrm = col*rm; 446807171c5SHong Zhang for (i=0; i<rm; i++) { /* over rows of R in those columns */ 447807171c5SHong Zhang r1 = r2 = r3 = r4 = 0.0; 448807171c5SHong Zhang n = r->i[i+1] - r->i[i]; 449807171c5SHong Zhang rj = r->j + r->i[i]; 450807171c5SHong Zhang ra = r->a + r->i[i]; 451807171c5SHong Zhang for (j=0; j<n; j++) { 452807171c5SHong Zhang r1 += (*ra)*c[*rj]; 453807171c5SHong Zhang r2 += (*ra)*c2[*rj]; 454807171c5SHong Zhang r3 += (*ra)*c3[*rj]; 455807171c5SHong Zhang r4 += (*ra++)*c4[*rj++]; 456807171c5SHong Zhang } 457807171c5SHong Zhang d[colrm + i] = r1; 458807171c5SHong Zhang d[colrm + rm + i] = r2; 459807171c5SHong Zhang d[colrm + rm2 + i] = r3; 460807171c5SHong Zhang d[colrm + rm3 + i] = r4; 461807171c5SHong Zhang } 462807171c5SHong Zhang } 463807171c5SHong Zhang for (; col<cn; col++) { /* over extra columns of C */ 464807171c5SHong Zhang for (i=0; i<am; i++) { /* over rows of A in those columns */ 465807171c5SHong Zhang r1 = 0.0; 466807171c5SHong Zhang n = a->i[i+1] - a->i[i]; 467807171c5SHong Zhang aj = a->j + a->i[i]; 468807171c5SHong Zhang aa = a->a + a->i[i]; 469807171c5SHong Zhang for (j=0; j<n; j++) { 470807171c5SHong Zhang r1 += (*aa++)*b1[*aj++]; 471807171c5SHong Zhang } 472807171c5SHong Zhang c[i] = r1; 473807171c5SHong Zhang } 474807171c5SHong Zhang b1 += bm; 475807171c5SHong Zhang 476807171c5SHong Zhang for (i=0; i<rm; i++) { /* over rows of R in those columns */ 477807171c5SHong Zhang r1 = 0.0; 478807171c5SHong Zhang n = r->i[i+1] - r->i[i]; 479807171c5SHong Zhang rj = r->j + r->i[i]; 480807171c5SHong Zhang ra = r->a + r->i[i]; 481807171c5SHong Zhang for (j=0; j<n; j++) { 482807171c5SHong Zhang r1 += (*ra++)*c[*rj++]; 483807171c5SHong Zhang } 484807171c5SHong Zhang d[col*rm + i] = r1; 485807171c5SHong Zhang } 486807171c5SHong Zhang } 487807171c5SHong Zhang ierr = PetscLogFlops(cn*2.0*(a->nz + r->nz));CHKERRQ(ierr); 488807171c5SHong Zhang 4898c778c55SBarry Smith ierr = MatDenseRestoreArray(B,&b);CHKERRQ(ierr); 4908c778c55SBarry Smith ierr = MatDenseRestoreArray(RAB,&d);CHKERRQ(ierr); 491807171c5SHong Zhang ierr = MatAssemblyBegin(RAB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 492807171c5SHong Zhang ierr = MatAssemblyEnd(RAB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 493807171c5SHong Zhang PetscFunctionReturn(0); 494807171c5SHong Zhang } 495807171c5SHong Zhang 496807171c5SHong Zhang #undef __FUNCT__ 497807171c5SHong Zhang #define __FUNCT__ "MatRARtNumeric_SeqAIJ_SeqAIJ" 498807171c5SHong Zhang PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ(Mat A,Mat R,Mat C) 499807171c5SHong Zhang { 500807171c5SHong Zhang PetscErrorCode ierr; 501807171c5SHong Zhang Mat_RARt *rart; 502807171c5SHong Zhang PetscContainer container; 503807171c5SHong Zhang MatTransposeColoring matcoloring; 5048d0a38e4SHong Zhang Mat Rt,RARt; 505807171c5SHong Zhang PetscLogDouble Mult_sp_den=0.0,app1=0.0,app2=0.0,t0,tf; 506807171c5SHong Zhang 507807171c5SHong Zhang PetscFunctionBegin; 508807171c5SHong Zhang ierr = PetscObjectQuery((PetscObject)C,"Mat_RARt",(PetscObject*)&container);CHKERRQ(ierr); 509807171c5SHong Zhang if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit"); 510807171c5SHong Zhang ierr = PetscContainerGetPointer(container,(void**)&rart);CHKERRQ(ierr); 511807171c5SHong Zhang 512807171c5SHong Zhang /* Get dense Rt by Apply MatTransposeColoring to R */ 513807171c5SHong Zhang matcoloring = rart->matcoloring; 514807171c5SHong Zhang Rt = rart->Rt; 5152205254eSKarl Rupp 5168563dfccSBarry Smith ierr = PetscTime(&t0);CHKERRQ(ierr); 517807171c5SHong Zhang ierr = MatTransColoringApplySpToDen(matcoloring,R,Rt);CHKERRQ(ierr); 5188563dfccSBarry Smith ierr = PetscTime(&tf);CHKERRQ(ierr); 519807171c5SHong Zhang app1 += tf - t0; 520807171c5SHong Zhang 521807171c5SHong Zhang /* Get dense RARt = R*A*Rt */ 5228563dfccSBarry Smith ierr = PetscTime(&t0);CHKERRQ(ierr); 523807171c5SHong Zhang RARt = rart->RARt; 524807171c5SHong Zhang ierr = MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense(R,A,Rt,RARt,rart->work);CHKERRQ(ierr); 5258563dfccSBarry Smith ierr = PetscTime(&tf);CHKERRQ(ierr); 5262205254eSKarl Rupp 527807171c5SHong Zhang Mult_sp_den += tf - t0; 528807171c5SHong Zhang 529807171c5SHong Zhang /* Recover C from C_dense */ 5308563dfccSBarry Smith ierr = PetscTime(&t0);CHKERRQ(ierr); 531807171c5SHong Zhang ierr = MatTransColoringApplyDenToSp(matcoloring,RARt,C);CHKERRQ(ierr); 5328563dfccSBarry Smith ierr = PetscTime(&tf);CHKERRQ(ierr); 5332205254eSKarl Rupp 534807171c5SHong Zhang app2 += tf - t0; 535807171c5SHong Zhang 536807171c5SHong Zhang #if defined(PETSC_USE_INFO) 537a2ea699eSBarry 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); 538807171c5SHong Zhang #endif 539807171c5SHong Zhang PetscFunctionReturn(0); 540807171c5SHong Zhang } 541807171c5SHong Zhang 542807171c5SHong Zhang #undef __FUNCT__ 543807171c5SHong Zhang #define __FUNCT__ "MatRARt_SeqAIJ_SeqAIJ" 544807171c5SHong Zhang PetscErrorCode MatRARt_SeqAIJ_SeqAIJ(Mat A,Mat R,MatReuse scall,PetscReal fill,Mat *C) 545807171c5SHong Zhang { 546807171c5SHong Zhang PetscErrorCode ierr; 547*22e94b5dSHong Zhang PetscBool usecoloring = PETSC_TRUE; 548807171c5SHong Zhang 549807171c5SHong Zhang PetscFunctionBegin; 550*22e94b5dSHong Zhang ierr = PetscOptionsGetBool(NULL,"-matrart_color",&usecoloring,NULL);CHKERRQ(ierr); 551*22e94b5dSHong Zhang 552*22e94b5dSHong Zhang if (!usecoloring) { 553*22e94b5dSHong Zhang Mat Rt; 554*22e94b5dSHong Zhang ierr = MatTranspose(R,MAT_INITIAL_MATRIX,&Rt);CHKERRQ(ierr); /* replace MAT_INITIAL_MATRIX with scall if !usecoloring is better */ 555*22e94b5dSHong Zhang ierr = MatMatMatMult(R,A,Rt,scall,fill,C);CHKERRQ(ierr); 556*22e94b5dSHong Zhang ierr = MatDestroy(&Rt);CHKERRQ(ierr); 557*22e94b5dSHong Zhang PetscFunctionReturn(0); 558*22e94b5dSHong Zhang } 559*22e94b5dSHong Zhang 560*22e94b5dSHong Zhang /* use coloring */ 561807171c5SHong Zhang if (scall == MAT_INITIAL_MATRIX) { 562807171c5SHong Zhang ierr = MatRARtSymbolic_SeqAIJ_SeqAIJ(A,R,fill,C);CHKERRQ(ierr); 563807171c5SHong Zhang } 564807171c5SHong Zhang ierr = MatRARtNumeric_SeqAIJ_SeqAIJ(A,R,*C);CHKERRQ(ierr); 565807171c5SHong Zhang PetscFunctionReturn(0); 566807171c5SHong Zhang } 567