1be1d678aSKris Buschelman 2d50806bdSBarry Smith /* 32499ec78SHong Zhang Defines matrix-matrix product routines for pairs of SeqAIJ matrices 4d50806bdSBarry Smith C = A * B 5d50806bdSBarry Smith */ 6d50806bdSBarry Smith 7c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/aij.h> /*I "petscmat.h" I*/ 8c6db04a5SJed Brown #include <../src/mat/utils/freespace.h> 9c6db04a5SJed Brown #include <petscbt.h> 10c6db04a5SJed Brown #include <../src/mat/impls/dense/seq/dense.h> /*I "petscmat.h" I*/ 116284ec50SHong Zhang 12ec01deb9SMatthew Knepley EXTERN_C_BEGIN 136284ec50SHong Zhang #undef __FUNCT__ 145c66b693SKris Buschelman #define __FUNCT__ "MatMatMult_SeqAIJ_SeqAIJ" 1538baddfdSBarry Smith PetscErrorCode MatMatMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 1638baddfdSBarry Smith { 17dfbe8321SBarry Smith PetscErrorCode ierr; 185c66b693SKris Buschelman 195c66b693SKris Buschelman PetscFunctionBegin; 2026be0446SHong Zhang if (scall == MAT_INITIAL_MATRIX){ 2153565b12SHong Zhang /* ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr); */ 2226be0446SHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr); 2353565b12SHong Zhang /* ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr); */ 2426be0446SHong Zhang } 2553565b12SHong Zhang /* ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr); */ 2601e47db4SHong Zhang ierr = (*(*C)->ops->matmultnumeric)(A,B,*C);CHKERRQ(ierr); 2753565b12SHong Zhang /* ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr); */ 285c66b693SKris Buschelman PetscFunctionReturn(0); 295c66b693SKris Buschelman } 30ec01deb9SMatthew Knepley EXTERN_C_END 311c24bd37SHong Zhang 3226be0446SHong Zhang #undef __FUNCT__ 33b561aa9dSHong Zhang #define __FUNCT__ "MatGetSymbolicMatMatMult_SeqAIJ_SeqAIJ" 34b561aa9dSHong Zhang PetscErrorCode MatGetSymbolicMatMatMult_SeqAIJ_SeqAIJ(PetscInt am,PetscInt *Ai,PetscInt *Aj,PetscInt bm,PetscInt bn,PetscInt *Bi,PetscInt *Bj,PetscReal fill,PetscInt *Ci[],PetscInt *Cj[],PetscInt *nspacedouble) 3558c24d83SHong Zhang { 36dfbe8321SBarry Smith PetscErrorCode ierr; 37a1a86e44SBarry Smith PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 38b561aa9dSHong Zhang PetscInt *ai=Ai,*aj=Aj,*bi=Bi,*bj=Bj,*bjj,*ci,*cj; 39b561aa9dSHong Zhang PetscInt i,j,anzi,brow,bnzj,cnzi,nlnk,*lnk,ndouble=0; 40be0fcf8dSHong Zhang PetscBT lnkbt; 4158c24d83SHong Zhang 4258c24d83SHong Zhang PetscFunctionBegin; 4358c24d83SHong Zhang /* Allocate ci array, arrays for fill computation and */ 4458c24d83SHong Zhang /* free space for accumulating nonzero column info */ 4538baddfdSBarry Smith ierr = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr); 4658c24d83SHong Zhang ci[0] = 0; 4758c24d83SHong Zhang 48be0fcf8dSHong Zhang /* create and initialize a linked list */ 49be0fcf8dSHong Zhang nlnk = bn+1; 50be0fcf8dSHong Zhang ierr = PetscLLCreate(bn,bn,nlnk,lnk,lnkbt);CHKERRQ(ierr); 5158c24d83SHong Zhang 52c5db241fSHong Zhang /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */ 53a1a86e44SBarry Smith ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr); 5458c24d83SHong Zhang current_space = free_space; 5558c24d83SHong Zhang 5658c24d83SHong Zhang /* Determine symbolic info for each row of the product: */ 5758c24d83SHong Zhang for (i=0;i<am;i++) { 5858c24d83SHong Zhang anzi = ai[i+1] - ai[i]; 5958c24d83SHong Zhang cnzi = 0; 602d09714cSHong Zhang j = anzi; 61b561aa9dSHong Zhang aj = Aj + ai[i]; 622d09714cSHong Zhang while (j){/* assume cols are almost in increasing order, starting from its end saves computation */ 632d09714cSHong Zhang j--; 64b561aa9dSHong Zhang brow = aj[j]; 6558c24d83SHong Zhang bnzj = bi[brow+1] - bi[brow]; 6658c24d83SHong Zhang bjj = bj + bi[brow]; 671c239cc6SHong Zhang /* add non-zero cols of B into the sorted linked list lnk */ 68be0fcf8dSHong Zhang ierr = PetscLLAdd(bnzj,bjj,bn,nlnk,lnk,lnkbt);CHKERRQ(ierr); 691c239cc6SHong Zhang cnzi += nlnk; 7058c24d83SHong Zhang } 7158c24d83SHong Zhang 7258c24d83SHong Zhang /* If free space is not available, make more free space */ 7358c24d83SHong Zhang /* Double the amount of total space in the list */ 7458c24d83SHong Zhang if (current_space->local_remaining<cnzi) { 754238b7adSHong Zhang ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 76b561aa9dSHong Zhang ndouble++; 7758c24d83SHong Zhang } 7858c24d83SHong Zhang 79c5db241fSHong Zhang /* Copy data into free space, then initialize lnk */ 80be0fcf8dSHong Zhang ierr = PetscLLClean(bn,bn,cnzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 81c5db241fSHong Zhang current_space->array += cnzi; 8258c24d83SHong Zhang current_space->local_used += cnzi; 8358c24d83SHong Zhang current_space->local_remaining -= cnzi; 8458c24d83SHong Zhang ci[i+1] = ci[i] + cnzi; 8558c24d83SHong Zhang } 8658c24d83SHong Zhang 8758c24d83SHong Zhang /* Column indices are in the list of free space */ 8858c24d83SHong Zhang /* Allocate space for cj, initialize cj, and */ 8958c24d83SHong Zhang /* destroy list of free space and other temporary array(s) */ 9038baddfdSBarry Smith ierr = PetscMalloc((ci[am]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr); 91a1a86e44SBarry Smith ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr); 92be0fcf8dSHong Zhang ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 9358c24d83SHong Zhang 94b561aa9dSHong Zhang *Ci = ci; 95b561aa9dSHong Zhang *Cj = cj; 96b561aa9dSHong Zhang *nspacedouble = ndouble; 97b561aa9dSHong Zhang PetscFunctionReturn(0); 98b561aa9dSHong Zhang } 99b561aa9dSHong Zhang 100b561aa9dSHong Zhang #undef __FUNCT__ 101b561aa9dSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ" 102b561aa9dSHong Zhang PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C) 103b561aa9dSHong Zhang { 104b561aa9dSHong Zhang PetscErrorCode ierr; 105b561aa9dSHong Zhang Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c; 106b561aa9dSHong Zhang PetscInt *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*ci,*cj; 107b561aa9dSHong Zhang PetscInt am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N,nspacedouble; 108b561aa9dSHong Zhang MatScalar *ca; 109b561aa9dSHong Zhang PetscReal afill; 11036ec6d2dSHong Zhang PetscBool dense_axpy; /* false: use sparse axpy; otherwise use dense axpy in MatMatMultNumeric_SeqAIJ_SeqAIJ() */ 111b561aa9dSHong Zhang 112b561aa9dSHong Zhang PetscFunctionBegin; 113b561aa9dSHong Zhang /* Get ci and cj */ 114b561aa9dSHong Zhang ierr = MatGetSymbolicMatMatMult_SeqAIJ_SeqAIJ(am,ai,aj,bm,bn,bi,bj,fill,&ci,&cj,&nspacedouble);CHKERRQ(ierr); 115b561aa9dSHong Zhang 11658c24d83SHong Zhang /* Allocate space for ca */ 11758c24d83SHong Zhang ierr = PetscMalloc((ci[am]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr); 11858c24d83SHong Zhang ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr); 11958c24d83SHong Zhang 12026be0446SHong Zhang /* put together the new symbolic matrix */ 1217adad957SLisandro Dalcin ierr = MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,am,bn,ci,cj,ca,C);CHKERRQ(ierr); 12258c24d83SHong Zhang 12358c24d83SHong Zhang /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 12458c24d83SHong Zhang /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */ 12558c24d83SHong Zhang c = (Mat_SeqAIJ *)((*C)->data); 126e6b907acSBarry Smith c->free_a = PETSC_TRUE; 127e6b907acSBarry Smith c->free_ij = PETSC_TRUE; 12858c24d83SHong Zhang c->nonew = 0; 1298cdbd757SHong Zhang (*C)->ops->matmult = MatMatMult_SeqAIJ_SeqAIJ; 13058c24d83SHong Zhang 13153565b12SHong Zhang /* Determine which MatMatMultNumeric_SeqAIJ_SeqAIJ() to be used */ 13236ec6d2dSHong Zhang dense_axpy = PETSC_TRUE; 13336ec6d2dSHong Zhang ierr = PetscOptionsGetBool(PETSC_NULL,"-matmatmult_denseaxpy",&dense_axpy,PETSC_NULL);CHKERRQ(ierr); 13436ec6d2dSHong Zhang if (dense_axpy){ 13536ec6d2dSHong Zhang ierr = PetscMalloc(bn*sizeof(PetscScalar),&c->matmult_abdense);CHKERRQ(ierr); 136c58ca2e3SHong Zhang ierr = PetscMemzero(c->matmult_abdense,dense_axpy*bn*sizeof(PetscScalar));CHKERRQ(ierr); 137c58ca2e3SHong Zhang (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ; /* fast, takes additional dense_axpy*bn*sizeof(PetscScalar) space */ 138c58ca2e3SHong Zhang } else { /* slower, but use less memory */ 139c58ca2e3SHong Zhang (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ_SparseAxpy; /* slower, less memory */ 140c58ca2e3SHong Zhang } 1410b7e3e3dSHong Zhang 1427212da7cSHong Zhang /* set MatInfo */ 1437212da7cSHong Zhang afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5; 144f2b054eeSHong Zhang if (afill < 1.0) afill = 1.0; 1457212da7cSHong Zhang c->maxnz = ci[am]; 1467212da7cSHong Zhang c->nz = ci[am]; 1477212da7cSHong Zhang (*C)->info.mallocs = nspacedouble; 1487212da7cSHong Zhang (*C)->info.fill_ratio_given = fill; 1497212da7cSHong Zhang (*C)->info.fill_ratio_needed = afill; 1507212da7cSHong Zhang 1517212da7cSHong Zhang #if defined(PETSC_USE_INFO) 1527212da7cSHong Zhang if (ci[am]) { 153f2b054eeSHong Zhang ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr); 154f2b054eeSHong Zhang ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr); 155f2b054eeSHong Zhang } else { 156f2b054eeSHong Zhang ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr); 157be0fcf8dSHong Zhang } 158f2b054eeSHong Zhang #endif 15958c24d83SHong Zhang PetscFunctionReturn(0); 16058c24d83SHong Zhang } 161d50806bdSBarry Smith 16226be0446SHong Zhang #undef __FUNCT__ 16326be0446SHong Zhang #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqAIJ" 164dfbe8321SBarry Smith PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C) 165d50806bdSBarry Smith { 166dfbe8321SBarry Smith PetscErrorCode ierr; 167d13dce4bSSatish Balay PetscLogDouble flops=0.0; 168d50806bdSBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 169d50806bdSBarry Smith Mat_SeqAIJ *b = (Mat_SeqAIJ *)B->data; 170d50806bdSBarry Smith Mat_SeqAIJ *c = (Mat_SeqAIJ *)C->data; 17138baddfdSBarry Smith PetscInt *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j; 172c58ca2e3SHong Zhang PetscInt am=A->rmap->n,cm=C->rmap->n; 173519eb980SHong Zhang PetscInt i,j,k,anzi,bnzi,cnzi,brow; 17436ec6d2dSHong Zhang PetscScalar *aa=a->a,*ba=b->a,*baj,*ca=c->a,valtmp; 1750b7e3e3dSHong Zhang PetscScalar *ab_dense=c->matmult_abdense; 176d50806bdSBarry Smith 177d50806bdSBarry Smith PetscFunctionBegin; 178c124e916SHong Zhang /* clean old values in C */ 179c124e916SHong Zhang ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr); 180d50806bdSBarry Smith /* Traverse A row-wise. */ 181d50806bdSBarry Smith /* Build the ith row in C by summing over nonzero columns in A, */ 182d50806bdSBarry Smith /* the rows of B corresponding to nonzeros of A. */ 183d50806bdSBarry Smith for (i=0; i<am; i++) { 184d50806bdSBarry Smith anzi = ai[i+1] - ai[i]; 185d50806bdSBarry Smith for (j=0; j<anzi; j++) { 186519eb980SHong Zhang brow = aj[j]; 187d50806bdSBarry Smith bnzi = bi[brow+1] - bi[brow]; 188d50806bdSBarry Smith bjj = bj + bi[brow]; 189d50806bdSBarry Smith baj = ba + bi[brow]; 190519eb980SHong Zhang /* perform dense axpy */ 19136ec6d2dSHong Zhang valtmp = aa[j]; 192519eb980SHong Zhang for (k=0; k<bnzi; k++) { 19336ec6d2dSHong Zhang ab_dense[bjj[k]] += valtmp*baj[k]; 194519eb980SHong Zhang } 195519eb980SHong Zhang flops += 2*bnzi; 196519eb980SHong Zhang } 197c58ca2e3SHong Zhang aj += anzi; aa += anzi; 198c58ca2e3SHong Zhang 199c58ca2e3SHong Zhang cnzi = ci[i+1] - ci[i]; 200519eb980SHong Zhang for (k=0; k<cnzi; k++) { 201519eb980SHong Zhang ca[k] += ab_dense[cj[k]]; 202519eb980SHong Zhang ab_dense[cj[k]] = 0.0; /* zero ab_dense */ 203519eb980SHong Zhang } 204519eb980SHong Zhang flops += cnzi; 205519eb980SHong Zhang cj += cnzi; ca += cnzi; 206519eb980SHong Zhang } 207c58ca2e3SHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 208c58ca2e3SHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 209c58ca2e3SHong Zhang ierr = PetscLogFlops(flops);CHKERRQ(ierr); 210c58ca2e3SHong Zhang C->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ; 211c58ca2e3SHong Zhang PetscFunctionReturn(0); 212c58ca2e3SHong Zhang } 213c58ca2e3SHong Zhang 214c58ca2e3SHong Zhang #undef __FUNCT__ 215c58ca2e3SHong Zhang #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqAIJ_SparseAxpy" 216c58ca2e3SHong Zhang PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_SparseAxpy(Mat A,Mat B,Mat C) 217c58ca2e3SHong Zhang { 218c58ca2e3SHong Zhang PetscErrorCode ierr; 219c58ca2e3SHong Zhang PetscLogDouble flops=0.0; 220c58ca2e3SHong Zhang Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 221c58ca2e3SHong Zhang Mat_SeqAIJ *b = (Mat_SeqAIJ *)B->data; 222c58ca2e3SHong Zhang Mat_SeqAIJ *c = (Mat_SeqAIJ *)C->data; 223c58ca2e3SHong Zhang PetscInt *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j; 224c58ca2e3SHong Zhang PetscInt am=A->rmap->N,cm=C->rmap->N; 225c58ca2e3SHong Zhang PetscInt i,j,k,anzi,bnzi,cnzi,brow; 22636ec6d2dSHong Zhang PetscScalar *aa=a->a,*ba=b->a,*baj,*ca=c->a,valtmp; 227c58ca2e3SHong Zhang PetscInt nextb; 228c58ca2e3SHong Zhang 229c58ca2e3SHong Zhang PetscFunctionBegin; 230c58ca2e3SHong Zhang /* clean old values in C */ 231c58ca2e3SHong Zhang ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr); 232c58ca2e3SHong Zhang /* Traverse A row-wise. */ 233c58ca2e3SHong Zhang /* Build the ith row in C by summing over nonzero columns in A, */ 234c58ca2e3SHong Zhang /* the rows of B corresponding to nonzeros of A. */ 235519eb980SHong Zhang for (i=0;i<am;i++) { 236519eb980SHong Zhang anzi = ai[i+1] - ai[i]; 237519eb980SHong Zhang cnzi = ci[i+1] - ci[i]; 238519eb980SHong Zhang for (j=0;j<anzi;j++) { 239519eb980SHong Zhang brow = aj[j]; 240519eb980SHong Zhang bnzi = bi[brow+1] - bi[brow]; 241519eb980SHong Zhang bjj = bj + bi[brow]; 242519eb980SHong Zhang baj = ba + bi[brow]; 243519eb980SHong Zhang /* perform sparse axpy */ 24436ec6d2dSHong Zhang valtmp = aa[j]; 245c124e916SHong Zhang nextb = 0; 246c124e916SHong Zhang for (k=0; nextb<bnzi; k++) { 247c124e916SHong Zhang if (cj[k] == bjj[nextb]){ /* ccol == bcol */ 24836ec6d2dSHong Zhang ca[k] += valtmp*baj[nextb++]; 249c124e916SHong Zhang } 250d50806bdSBarry Smith } 251d50806bdSBarry Smith flops += 2*bnzi; 252d50806bdSBarry Smith } 253519eb980SHong Zhang aj += anzi; aa += anzi; 254519eb980SHong Zhang cj += cnzi; ca += cnzi; 255d50806bdSBarry Smith } 256c58ca2e3SHong Zhang 257716bacf3SKris Buschelman ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 258716bacf3SKris Buschelman ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 259d50806bdSBarry Smith ierr = PetscLogFlops(flops);CHKERRQ(ierr); 260d50806bdSBarry Smith PetscFunctionReturn(0); 261d50806bdSBarry Smith } 262bc011b1eSHong Zhang 263d2853540SHong Zhang /* This routine is not used. Should be removed! */ 264bc011b1eSHong Zhang #undef __FUNCT__ 2656fc122caSHong Zhang #define __FUNCT__ "MatMatTransposeMult_SeqAIJ_SeqAIJ" 2666fc122caSHong Zhang PetscErrorCode MatMatTransposeMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 2675df89d91SHong Zhang { 268bc011b1eSHong Zhang PetscErrorCode ierr; 269bc011b1eSHong Zhang 270bc011b1eSHong Zhang PetscFunctionBegin; 271bc011b1eSHong Zhang if (scall == MAT_INITIAL_MATRIX){ 2726fc122caSHong Zhang ierr = MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr); 273bc011b1eSHong Zhang } 2746fc122caSHong Zhang ierr = MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr); 275bc011b1eSHong Zhang PetscFunctionReturn(0); 276bc011b1eSHong Zhang } 277bc011b1eSHong Zhang 278bc011b1eSHong Zhang #undef __FUNCT__ 279e286af10SHong Zhang #define __FUNCT__ "PetscContainerDestroy_Mat_MatMatTransMult" 280e286af10SHong Zhang PetscErrorCode PetscContainerDestroy_Mat_MatMatTransMult(void *ptr) 2812128a86cSHong Zhang { 2822128a86cSHong Zhang PetscErrorCode ierr; 283e286af10SHong Zhang Mat_MatMatTransMult *multtrans=(Mat_MatMatTransMult*)ptr; 2842128a86cSHong Zhang 2852128a86cSHong Zhang PetscFunctionBegin; 286b9af6bddSHong Zhang ierr = MatTransposeColoringDestroy(&multtrans->matcoloring);CHKERRQ(ierr); 2872128a86cSHong Zhang ierr = MatDestroy(&multtrans->Bt_den);CHKERRQ(ierr); 2882128a86cSHong Zhang ierr = MatDestroy(&multtrans->ABt_den);CHKERRQ(ierr); 2892128a86cSHong Zhang ierr = PetscFree(multtrans);CHKERRQ(ierr); 2902128a86cSHong Zhang PetscFunctionReturn(0); 2912128a86cSHong Zhang } 2922128a86cSHong Zhang 2932128a86cSHong Zhang #undef __FUNCT__ 2942128a86cSHong Zhang #define __FUNCT__ "MatDestroy_SeqAIJ_MatMatMultTrans" 2952128a86cSHong Zhang PetscErrorCode MatDestroy_SeqAIJ_MatMatMultTrans(Mat A) 2962128a86cSHong Zhang { 2972128a86cSHong Zhang PetscErrorCode ierr; 2982128a86cSHong Zhang PetscContainer container; 299e286af10SHong Zhang Mat_MatMatTransMult *multtrans=PETSC_NULL; 3002128a86cSHong Zhang 3012128a86cSHong Zhang PetscFunctionBegin; 302e286af10SHong Zhang ierr = PetscObjectQuery((PetscObject)A,"Mat_MatMatTransMult",(PetscObject *)&container);CHKERRQ(ierr); 3032128a86cSHong Zhang if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit"); 3042128a86cSHong Zhang ierr = PetscContainerGetPointer(container,(void **)&multtrans);CHKERRQ(ierr); 3052128a86cSHong Zhang A->ops->destroy = multtrans->destroy; 3062128a86cSHong Zhang if (A->ops->destroy) { 3072128a86cSHong Zhang ierr = (*A->ops->destroy)(A);CHKERRQ(ierr); 3082128a86cSHong Zhang } 309e286af10SHong Zhang ierr = PetscObjectCompose((PetscObject)A,"Mat_MatMatTransMult",0);CHKERRQ(ierr); 3102128a86cSHong Zhang PetscFunctionReturn(0); 3112128a86cSHong Zhang } 3122128a86cSHong Zhang 3132128a86cSHong Zhang #undef __FUNCT__ 3146fc122caSHong Zhang #define __FUNCT__ "MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ" 3156fc122caSHong Zhang PetscErrorCode MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C) 316bc011b1eSHong Zhang { 3175df89d91SHong Zhang PetscErrorCode ierr; 31881d82fe4SHong Zhang Mat Bt; 31981d82fe4SHong Zhang PetscInt *bti,*btj; 320e286af10SHong Zhang Mat_MatMatTransMult *multtrans; 3212128a86cSHong Zhang PetscContainer container; 322d2853540SHong Zhang PetscLogDouble t0,tf,etime2=0.0; 323d2853540SHong Zhang 32481d82fe4SHong Zhang PetscFunctionBegin; 325d2853540SHong Zhang ierr = PetscGetTime(&t0);CHKERRQ(ierr); 32681d82fe4SHong Zhang /* create symbolic Bt */ 32781d82fe4SHong Zhang ierr = MatGetSymbolicTranspose_SeqAIJ(B,&bti,&btj);CHKERRQ(ierr); 32881d82fe4SHong Zhang ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,B->cmap->n,B->rmap->n,bti,btj,PETSC_NULL,&Bt);CHKERRQ(ierr); 32981d82fe4SHong Zhang 33081d82fe4SHong Zhang /* get symbolic C=A*Bt */ 33181d82fe4SHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,Bt,fill,C);CHKERRQ(ierr); 33281d82fe4SHong Zhang 3332128a86cSHong Zhang /* create a supporting struct for reuse intermidiate dense matrices with matcoloring */ 334e286af10SHong Zhang ierr = PetscNew(Mat_MatMatTransMult,&multtrans);CHKERRQ(ierr); 3352128a86cSHong Zhang 3362128a86cSHong Zhang /* attach the supporting struct to C */ 3372128a86cSHong Zhang ierr = PetscContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr); 3382128a86cSHong Zhang ierr = PetscContainerSetPointer(container,multtrans);CHKERRQ(ierr); 339e286af10SHong Zhang ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_Mat_MatMatTransMult);CHKERRQ(ierr); 340e286af10SHong Zhang ierr = PetscObjectCompose((PetscObject)(*C),"Mat_MatMatTransMult",(PetscObject)container);CHKERRQ(ierr); 3412128a86cSHong Zhang ierr = PetscContainerDestroy(&container);CHKERRQ(ierr); 3422128a86cSHong Zhang 3432128a86cSHong Zhang multtrans->usecoloring = PETSC_FALSE; 3442128a86cSHong Zhang multtrans->destroy = (*C)->ops->destroy; 3452128a86cSHong Zhang (*C)->ops->destroy = MatDestroy_SeqAIJ_MatMatMultTrans; 3462128a86cSHong Zhang 347d2853540SHong Zhang ierr = PetscGetTime(&tf);CHKERRQ(ierr); 348d2853540SHong Zhang etime2 += tf - t0; 349d2853540SHong Zhang 350f400d4cbSHong Zhang ierr = PetscOptionsGetBool(PETSC_NULL,"-matmattransmult_color",&multtrans->usecoloring,PETSC_NULL);CHKERRQ(ierr); 3512128a86cSHong Zhang if (multtrans->usecoloring){ 352b9af6bddSHong Zhang /* Create MatTransposeColoring from symbolic C=A*B^T */ 353b9af6bddSHong Zhang MatTransposeColoring matcoloring; 3542128a86cSHong Zhang ISColoring iscoloring; 3552128a86cSHong Zhang Mat Bt_dense,C_dense; 356d2853540SHong Zhang PetscLogDouble etime0=0.0,etime01=0.0,etime1=0.0; 357e8354b3eSHong Zhang 358e8354b3eSHong Zhang ierr = PetscGetTime(&t0);CHKERRQ(ierr); 359502de53fSHong Zhang ierr = MatGetColoring(*C,MATCOLORINGLF,&iscoloring);CHKERRQ(ierr); 360d2853540SHong Zhang ierr = PetscGetTime(&tf);CHKERRQ(ierr); 361d2853540SHong Zhang etime0 += tf - t0; 362d2853540SHong Zhang 363d2853540SHong Zhang ierr = PetscGetTime(&t0);CHKERRQ(ierr); 364b9af6bddSHong Zhang ierr = MatTransposeColoringCreate(*C,iscoloring,&matcoloring);CHKERRQ(ierr); 3652128a86cSHong Zhang multtrans->matcoloring = matcoloring; 3662128a86cSHong Zhang ierr = ISColoringDestroy(&iscoloring);CHKERRQ(ierr); 367e8354b3eSHong Zhang ierr = PetscGetTime(&tf);CHKERRQ(ierr); 368d2853540SHong Zhang etime01 += tf - t0; 3692128a86cSHong Zhang 370e8354b3eSHong Zhang ierr = PetscGetTime(&t0);CHKERRQ(ierr); 3712128a86cSHong Zhang /* Create Bt_dense and C_dense = A*Bt_dense */ 3722128a86cSHong Zhang ierr = MatCreate(PETSC_COMM_SELF,&Bt_dense);CHKERRQ(ierr); 3732128a86cSHong Zhang ierr = MatSetSizes(Bt_dense,A->cmap->n,matcoloring->ncolors,A->cmap->n,matcoloring->ncolors);CHKERRQ(ierr); 3742128a86cSHong Zhang ierr = MatSetType(Bt_dense,MATSEQDENSE);CHKERRQ(ierr); 3752128a86cSHong Zhang ierr = MatSeqDenseSetPreallocation(Bt_dense,PETSC_NULL);CHKERRQ(ierr); 3762128a86cSHong Zhang Bt_dense->assembled = PETSC_TRUE; 3772128a86cSHong Zhang multtrans->Bt_den = Bt_dense; 3782128a86cSHong Zhang 3792128a86cSHong Zhang ierr = MatCreate(PETSC_COMM_SELF,&C_dense);CHKERRQ(ierr); 3802128a86cSHong Zhang ierr = MatSetSizes(C_dense,A->rmap->n,matcoloring->ncolors,A->rmap->n,matcoloring->ncolors);CHKERRQ(ierr); 3812128a86cSHong Zhang ierr = MatSetType(C_dense,MATSEQDENSE);CHKERRQ(ierr); 3822128a86cSHong Zhang ierr = MatSeqDenseSetPreallocation(C_dense,PETSC_NULL);CHKERRQ(ierr); 3832128a86cSHong Zhang Bt_dense->assembled = PETSC_TRUE; 3842128a86cSHong Zhang multtrans->ABt_den = C_dense; 385e8354b3eSHong Zhang ierr = PetscGetTime(&tf);CHKERRQ(ierr); 386e8354b3eSHong Zhang etime1 += tf - t0; 387f94ccd6cSHong Zhang 388f94ccd6cSHong Zhang #if defined(PETSC_USE_INFO) 3891ce71dffSSatish Balay { 390f94ccd6cSHong Zhang Mat_SeqAIJ *c=(Mat_SeqAIJ*)(*C)->data; 391f94ccd6cSHong Zhang ierr = PetscInfo5(*C,"Bt_dense: %D,%D; Cnz %D / (cm*ncolors %D) = %g\n",A->cmap->n,matcoloring->ncolors,c->nz,A->rmap->n*matcoloring->ncolors,(PetscReal)(c->nz)/(A->rmap->n*matcoloring->ncolors)); 392f94ccd6cSHong Zhang ierr = PetscInfo5(*C,"Sym = GetColor %g + ColorCreate %g + MatDenseCreate %g + non-colorSym %g = %g\n",etime0,etime01,etime1,etime2,etime0+etime01+etime1+etime2); 3931ce71dffSSatish Balay } 394f94ccd6cSHong Zhang #endif 3952128a86cSHong Zhang } 396*e99be685SHong Zhang /* clean up */ 397*e99be685SHong Zhang ierr = MatDestroy(&Bt);CHKERRQ(ierr); 398*e99be685SHong Zhang ierr = MatRestoreSymbolicTranspose_SeqAIJ(B,&bti,&btj);CHKERRQ(ierr); 3992128a86cSHong Zhang 400f94ccd6cSHong Zhang 401f94ccd6cSHong Zhang 40281d82fe4SHong Zhang #if defined(INEFFICIENT_ALGORITHM) 40381d82fe4SHong Zhang /* The algorithm below computes am*bm sparse inner-product - inefficient! It will be deleted later. */ 4041fa1209cSHong Zhang PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 4051fa1209cSHong Zhang Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c; 4061fa1209cSHong Zhang PetscInt *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*ci,*cj,*acol,*bcol; 4071fa1209cSHong Zhang PetscInt am=A->rmap->N,bm=B->rmap->N; 4081fa1209cSHong Zhang PetscInt i,j,anzi,bnzj,cnzi,nlnk,*lnk,nspacedouble=0,ka,kb,index[1]; 4091fa1209cSHong Zhang MatScalar *ca; 4101fa1209cSHong Zhang PetscBT lnkbt; 4111fa1209cSHong Zhang PetscReal afill; 4125df89d91SHong Zhang 4131fa1209cSHong Zhang /* Allocate row pointer array ci */ 4141fa1209cSHong Zhang ierr = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr); 4151fa1209cSHong Zhang ci[0] = 0; 4161fa1209cSHong Zhang 4171fa1209cSHong Zhang /* Create and initialize a linked list for C columns */ 4181fa1209cSHong Zhang nlnk = bm+1; 4191fa1209cSHong Zhang ierr = PetscLLCreate(bm,bm,nlnk,lnk,lnkbt);CHKERRQ(ierr); 4201fa1209cSHong Zhang 4211fa1209cSHong Zhang /* Initial FreeSpace with size fill*(nnz(A)+nnz(B)) */ 4221fa1209cSHong Zhang ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr); 4231fa1209cSHong Zhang current_space = free_space; 4241fa1209cSHong Zhang 4251fa1209cSHong Zhang /* Determine symbolic info for each row of the product A*B^T: */ 4261fa1209cSHong Zhang for (i=0; i<am; i++) { 4271fa1209cSHong Zhang anzi = ai[i+1] - ai[i]; 4281fa1209cSHong Zhang cnzi = 0; 4291fa1209cSHong Zhang acol = aj + ai[i]; 4301fa1209cSHong Zhang for (j=0; j<bm; j++){ 4311fa1209cSHong Zhang bnzj = bi[j+1] - bi[j]; 4321fa1209cSHong Zhang bcol= bj + bi[j]; 43381d82fe4SHong Zhang /* sparse inner-product c(i,j)=A[i,:]*B[j,:]^T */ 4341fa1209cSHong Zhang ka = 0; kb = 0; 4351fa1209cSHong Zhang while (ka < anzi && kb < bnzj){ 43681d82fe4SHong Zhang while (acol[ka] < bcol[kb] && ka < anzi) ka++; 43781d82fe4SHong Zhang if (ka == anzi) break; 43881d82fe4SHong Zhang while (acol[ka] > bcol[kb] && kb < bnzj) kb++; 43981d82fe4SHong Zhang if (kb == bnzj) break; 4401fa1209cSHong Zhang if (acol[ka] == bcol[kb]){ /* add nonzero c(i,j) to lnk */ 4411fa1209cSHong Zhang index[0] = j; 4421fa1209cSHong Zhang ierr = PetscLLAdd(1,index,bm,nlnk,lnk,lnkbt);CHKERRQ(ierr); 4431fa1209cSHong Zhang cnzi++; 4441fa1209cSHong Zhang break; 4451fa1209cSHong Zhang } 4461fa1209cSHong Zhang } 4471fa1209cSHong Zhang } 4481fa1209cSHong Zhang 4491fa1209cSHong Zhang /* If free space is not available, make more free space */ 4501fa1209cSHong Zhang /* Double the amount of total space in the list */ 4511fa1209cSHong Zhang if (current_space->local_remaining<cnzi) { 4521fa1209cSHong Zhang ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 4531fa1209cSHong Zhang nspacedouble++; 4541fa1209cSHong Zhang } 4551fa1209cSHong Zhang 4561fa1209cSHong Zhang /* Copy data into free space, then initialize lnk */ 4571fa1209cSHong Zhang ierr = PetscLLClean(bm,bm,cnzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 4581fa1209cSHong Zhang current_space->array += cnzi; 4591fa1209cSHong Zhang current_space->local_used += cnzi; 4601fa1209cSHong Zhang current_space->local_remaining -= cnzi; 4611fa1209cSHong Zhang 4621fa1209cSHong Zhang ci[i+1] = ci[i] + cnzi; 4631fa1209cSHong Zhang } 4641fa1209cSHong Zhang 4651fa1209cSHong Zhang 4661fa1209cSHong Zhang /* Column indices are in the list of free space. 4671fa1209cSHong Zhang Allocate array cj, copy column indices to cj, and destroy list of free space */ 4681fa1209cSHong Zhang ierr = PetscMalloc((ci[am]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr); 4691fa1209cSHong Zhang ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr); 4701fa1209cSHong Zhang ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 4711fa1209cSHong Zhang 4721fa1209cSHong Zhang /* Allocate space for ca */ 4731fa1209cSHong Zhang ierr = PetscMalloc((ci[am]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr); 4741fa1209cSHong Zhang ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr); 4751fa1209cSHong Zhang 4761fa1209cSHong Zhang /* put together the new symbolic matrix */ 4771fa1209cSHong Zhang ierr = MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,am,bm,ci,cj,ca,C);CHKERRQ(ierr); 4781fa1209cSHong Zhang 4791fa1209cSHong Zhang /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 4801fa1209cSHong Zhang /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */ 4811fa1209cSHong Zhang c = (Mat_SeqAIJ *)((*C)->data); 4821fa1209cSHong Zhang c->free_a = PETSC_TRUE; 4831fa1209cSHong Zhang c->free_ij = PETSC_TRUE; 4841fa1209cSHong Zhang c->nonew = 0; 4851fa1209cSHong Zhang 4861fa1209cSHong Zhang /* set MatInfo */ 4871fa1209cSHong Zhang afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5; 4881fa1209cSHong Zhang if (afill < 1.0) afill = 1.0; 4891fa1209cSHong Zhang c->maxnz = ci[am]; 4901fa1209cSHong Zhang c->nz = ci[am]; 4911fa1209cSHong Zhang (*C)->info.mallocs = nspacedouble; 4921fa1209cSHong Zhang (*C)->info.fill_ratio_given = fill; 4931fa1209cSHong Zhang (*C)->info.fill_ratio_needed = afill; 4941fa1209cSHong Zhang 4951fa1209cSHong Zhang #if defined(PETSC_USE_INFO) 4961fa1209cSHong Zhang if (ci[am]) { 4971fa1209cSHong Zhang ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr); 4986fc122caSHong Zhang ierr = PetscInfo1((*C),"Use MatMatTransposeMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr); 4991fa1209cSHong Zhang } else { 5001fa1209cSHong Zhang ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr); 5011fa1209cSHong Zhang } 5021fa1209cSHong Zhang #endif 50381d82fe4SHong Zhang #endif 5045df89d91SHong Zhang PetscFunctionReturn(0); 5055df89d91SHong Zhang } 5065df89d91SHong Zhang 5076973769bSHong Zhang /* #define USE_ARRAY - for sparse dot product. Slower than !USE_ARRAY */ 5085df89d91SHong Zhang #undef __FUNCT__ 5096fc122caSHong Zhang #define __FUNCT__ "MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ" 5106fc122caSHong Zhang PetscErrorCode MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C) 5115df89d91SHong Zhang { 5125df89d91SHong Zhang PetscErrorCode ierr; 5135df89d91SHong Zhang Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data; 514e2cac8caSJed Brown PetscInt *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,anzi,bnzj,nexta,nextb,*acol,*bcol,brow; 51581d82fe4SHong Zhang PetscInt cm=C->rmap->n,*ci=c->i,*cj=c->j,i,j,cnzi,*ccol; 5165df89d91SHong Zhang PetscLogDouble flops=0.0; 5176973769bSHong Zhang MatScalar *aa=a->a,*aval,*ba=b->a,*bval,*ca=c->a,*cval; 518e286af10SHong Zhang Mat_MatMatTransMult *multtrans; 5192128a86cSHong Zhang PetscContainer container; 5206973769bSHong Zhang #if defined(USE_ARRAY) 5216973769bSHong Zhang MatScalar *spdot; 5226973769bSHong Zhang #endif 5235df89d91SHong Zhang 5245df89d91SHong Zhang PetscFunctionBegin; 525e286af10SHong Zhang ierr = PetscObjectQuery((PetscObject)C,"Mat_MatMatTransMult",(PetscObject *)&container);CHKERRQ(ierr); 5262128a86cSHong Zhang if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit"); 5272128a86cSHong Zhang ierr = PetscContainerGetPointer(container,(void **)&multtrans);CHKERRQ(ierr); 5282128a86cSHong Zhang if (multtrans->usecoloring){ 529b9af6bddSHong Zhang MatTransposeColoring matcoloring = multtrans->matcoloring; 530c8db22f6SHong Zhang Mat Bt_dense; 531c8db22f6SHong Zhang PetscInt m,n; 532b2d2b04fSHong Zhang PetscLogDouble t0,tf,etime0=0.0,etime1=0.0,etime2=0.0; 533a0b95ffbSSatish Balay Mat C_dense = multtrans->ABt_den; 534a0b95ffbSSatish Balay 5352128a86cSHong Zhang Bt_dense = multtrans->Bt_den; 536c8db22f6SHong Zhang ierr = MatGetLocalSize(Bt_dense,&m,&n);CHKERRQ(ierr); 537c8db22f6SHong Zhang 538b9af6bddSHong Zhang /* Get Bt_dense by Apply MatTransposeColoring to B */ 5392128a86cSHong Zhang ierr = PetscGetTime(&t0);CHKERRQ(ierr); 540b9af6bddSHong Zhang ierr = MatTransColoringApplySpToDen(matcoloring,B,Bt_dense);CHKERRQ(ierr); 5412128a86cSHong Zhang ierr = PetscGetTime(&tf);CHKERRQ(ierr); 542b2d2b04fSHong Zhang etime0 += tf - t0; 543c8db22f6SHong Zhang 544c8db22f6SHong Zhang /* C_dense = A*Bt_dense */ 5452128a86cSHong Zhang ierr = PetscGetTime(&t0);CHKERRQ(ierr); 5462128a86cSHong Zhang ierr = MatMatMultNumeric_SeqAIJ_SeqDense(A,Bt_dense,C_dense);CHKERRQ(ierr); 5472128a86cSHong Zhang ierr = PetscGetTime(&tf);CHKERRQ(ierr); 5482128a86cSHong Zhang etime2 += tf - t0; 549c8db22f6SHong Zhang 5502128a86cSHong Zhang /* Recover C from C_dense */ 5512128a86cSHong Zhang ierr = PetscGetTime(&t0);CHKERRQ(ierr); 552b9af6bddSHong Zhang ierr = MatTransColoringApplyDenToSp(matcoloring,C_dense,C);CHKERRQ(ierr); 5532128a86cSHong Zhang ierr = PetscGetTime(&tf);CHKERRQ(ierr); 5542128a86cSHong Zhang etime1 += tf - t0; 555f94ccd6cSHong Zhang #if defined(PETSC_USE_INFO) 556f94ccd6cSHong Zhang ierr = PetscInfo4(C,"Num = ColoringApply: %g %g + Mult_sp_dense %g = %g\n",etime0,etime1,etime2,etime0+etime1+etime2); 557f94ccd6cSHong Zhang #endif 558c8db22f6SHong Zhang PetscFunctionReturn(0); 559c8db22f6SHong Zhang } 560c8db22f6SHong Zhang 5616973769bSHong Zhang #if defined(USE_ARRAY) 5626973769bSHong Zhang /* allocate an array for implementing sparse inner-product */ 563e2cac8caSJed Brown ierr = PetscMalloc((A->cmap->n+1)*sizeof(MatScalar),&spdot);CHKERRQ(ierr); 564e2cac8caSJed Brown ierr = PetscMemzero(spdot,(A->cmap->n+1)*sizeof(MatScalar));CHKERRQ(ierr); 5656973769bSHong Zhang #endif 5666973769bSHong Zhang 56781d82fe4SHong Zhang /* clear old values in C */ 56881d82fe4SHong Zhang ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr); 5695df89d91SHong Zhang 5701fa1209cSHong Zhang for (i=0; i<cm; i++) { 57181d82fe4SHong Zhang anzi = ai[i+1] - ai[i]; 57281d82fe4SHong Zhang acol = aj + ai[i]; 5736973769bSHong Zhang aval = aa + ai[i]; 5741fa1209cSHong Zhang cnzi = ci[i+1] - ci[i]; 5751fa1209cSHong Zhang ccol = cj + ci[i]; 5766973769bSHong Zhang cval = ca + ci[i]; 5771fa1209cSHong Zhang for (j=0; j<cnzi; j++){ 57881d82fe4SHong Zhang brow = ccol[j]; 57981d82fe4SHong Zhang bnzj = bi[brow+1] - bi[brow]; 58081d82fe4SHong Zhang bcol = bj + bi[brow]; 5816973769bSHong Zhang bval = ba + bi[brow]; 5826973769bSHong Zhang 58381d82fe4SHong Zhang /* perform sparse inner-product c(i,j)=A[i,:]*B[j,:]^T */ 5846973769bSHong Zhang #if defined(USE_ARRAY) 5856973769bSHong Zhang /* put ba to spdot array */ 5866973769bSHong Zhang for (nextb=0; nextb<bnzj; nextb++) spdot[bcol[nextb]] = bval[nextb]; 5876973769bSHong Zhang /* c(i,j)=A[i,:]*B[j,:]^T */ 5886973769bSHong Zhang for (nexta=0; nexta<anzi; nexta++){ 5896973769bSHong Zhang cval[j] += spdot[acol[nexta]]*aval[nexta]; 5906973769bSHong Zhang } 5916973769bSHong Zhang /* zero spdot array */ 5926973769bSHong Zhang for (nextb=0; nextb<bnzj; nextb++) spdot[bcol[nextb]] = 0.0; 5936973769bSHong Zhang #else 59481d82fe4SHong Zhang nexta = 0; nextb = 0; 59581d82fe4SHong Zhang while (nexta<anzi && nextb<bnzj){ 59681d82fe4SHong Zhang while (acol[nexta] < bcol[nextb] && nexta < anzi) nexta++; 59781d82fe4SHong Zhang if (nexta == anzi) break; 59881d82fe4SHong Zhang while (acol[nexta] > bcol[nextb] && nextb < bnzj) nextb++; 59981d82fe4SHong Zhang if (nextb == bnzj) break; 60081d82fe4SHong Zhang if (acol[nexta] == bcol[nextb]){ 6016973769bSHong Zhang cval[j] += aval[nexta]*bval[nextb]; 60281d82fe4SHong Zhang nexta++; nextb++; 60381d82fe4SHong Zhang flops += 2; 6041fa1209cSHong Zhang } 6051fa1209cSHong Zhang } 6066973769bSHong Zhang #endif 60781d82fe4SHong Zhang } 60881d82fe4SHong Zhang } 60981d82fe4SHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 61081d82fe4SHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 61181d82fe4SHong Zhang ierr = PetscLogFlops(flops);CHKERRQ(ierr); 6126973769bSHong Zhang #if defined(USE_ARRAY) 6136973769bSHong Zhang ierr = PetscFree(spdot); 6146973769bSHong Zhang #endif 6155df89d91SHong Zhang PetscFunctionReturn(0); 6165df89d91SHong Zhang } 6175df89d91SHong Zhang 6185df89d91SHong Zhang #undef __FUNCT__ 61975648e8dSHong Zhang #define __FUNCT__ "MatTransposeMatMult_SeqAIJ_SeqAIJ" 62075648e8dSHong Zhang PetscErrorCode MatTransposeMatMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) { 6215df89d91SHong Zhang PetscErrorCode ierr; 6225df89d91SHong Zhang 6235df89d91SHong Zhang PetscFunctionBegin; 6245df89d91SHong Zhang if (scall == MAT_INITIAL_MATRIX){ 62575648e8dSHong Zhang ierr = MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr); 6265df89d91SHong Zhang } 62775648e8dSHong Zhang ierr = MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr); 6285df89d91SHong Zhang PetscFunctionReturn(0); 6295df89d91SHong Zhang } 6305df89d91SHong Zhang 6315df89d91SHong Zhang #undef __FUNCT__ 63275648e8dSHong Zhang #define __FUNCT__ "MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ" 63375648e8dSHong Zhang PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C) 6345df89d91SHong Zhang { 635bc011b1eSHong Zhang PetscErrorCode ierr; 636bc011b1eSHong Zhang Mat At; 63738baddfdSBarry Smith PetscInt *ati,*atj; 638bc011b1eSHong Zhang 639bc011b1eSHong Zhang PetscFunctionBegin; 640bc011b1eSHong Zhang /* create symbolic At */ 641bc011b1eSHong Zhang ierr = MatGetSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr); 642d0f46423SBarry Smith ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,A->cmap->n,A->rmap->n,ati,atj,PETSC_NULL,&At);CHKERRQ(ierr); 643bc011b1eSHong Zhang 644bc011b1eSHong Zhang /* get symbolic C=At*B */ 645bc011b1eSHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(At,B,fill,C);CHKERRQ(ierr); 646bc011b1eSHong Zhang 647bc011b1eSHong Zhang /* clean up */ 6486bf464f9SBarry Smith ierr = MatDestroy(&At);CHKERRQ(ierr); 649bc011b1eSHong Zhang ierr = MatRestoreSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr); 650bc011b1eSHong Zhang PetscFunctionReturn(0); 651bc011b1eSHong Zhang } 652bc011b1eSHong Zhang 653bc011b1eSHong Zhang #undef __FUNCT__ 65475648e8dSHong Zhang #define __FUNCT__ "MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ" 65575648e8dSHong Zhang PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C) 656bc011b1eSHong Zhang { 657bc011b1eSHong Zhang PetscErrorCode ierr; 6580fbc74f4SHong Zhang Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data; 659d0f46423SBarry Smith PetscInt am=A->rmap->n,anzi,*ai=a->i,*aj=a->j,*bi=b->i,*bj,bnzi,nextb; 660d0f46423SBarry Smith PetscInt cm=C->rmap->n,*ci=c->i,*cj=c->j,crow,*cjj,i,j,k; 661d13dce4bSSatish Balay PetscLogDouble flops=0.0; 6620fbc74f4SHong Zhang MatScalar *aa=a->a,*ba,*ca=c->a,*caj; 663bc011b1eSHong Zhang 664bc011b1eSHong Zhang PetscFunctionBegin; 665bc011b1eSHong Zhang /* clear old values in C */ 666bc011b1eSHong Zhang ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr); 667bc011b1eSHong Zhang 668bc011b1eSHong Zhang /* compute A^T*B using outer product (A^T)[:,i]*B[i,:] */ 669bc011b1eSHong Zhang for (i=0;i<am;i++) { 670bc011b1eSHong Zhang bj = b->j + bi[i]; 671bc011b1eSHong Zhang ba = b->a + bi[i]; 672bc011b1eSHong Zhang bnzi = bi[i+1] - bi[i]; 673bc011b1eSHong Zhang anzi = ai[i+1] - ai[i]; 674bc011b1eSHong Zhang for (j=0; j<anzi; j++) { 675bc011b1eSHong Zhang nextb = 0; 6760fbc74f4SHong Zhang crow = *aj++; 677bc011b1eSHong Zhang cjj = cj + ci[crow]; 678bc011b1eSHong Zhang caj = ca + ci[crow]; 679bc011b1eSHong Zhang /* perform sparse axpy operation. Note cjj includes bj. */ 680bc011b1eSHong Zhang for (k=0; nextb<bnzi; k++) { 6810fbc74f4SHong Zhang if (cjj[k] == *(bj+nextb)) { /* ccol == bcol */ 6820fbc74f4SHong Zhang caj[k] += (*aa)*(*(ba+nextb)); 683bc011b1eSHong Zhang nextb++; 684bc011b1eSHong Zhang } 685bc011b1eSHong Zhang } 686bc011b1eSHong Zhang flops += 2*bnzi; 6870fbc74f4SHong Zhang aa++; 688bc011b1eSHong Zhang } 689bc011b1eSHong Zhang } 690bc011b1eSHong Zhang 691bc011b1eSHong Zhang /* Assemble the final matrix and clean up */ 692bc011b1eSHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 693bc011b1eSHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 694bc011b1eSHong Zhang ierr = PetscLogFlops(flops);CHKERRQ(ierr); 695bc011b1eSHong Zhang PetscFunctionReturn(0); 696bc011b1eSHong Zhang } 6979123193aSHong Zhang 698ec01deb9SMatthew Knepley EXTERN_C_BEGIN 6999123193aSHong Zhang #undef __FUNCT__ 7009123193aSHong Zhang #define __FUNCT__ "MatMatMult_SeqAIJ_SeqDense" 7019123193aSHong Zhang PetscErrorCode MatMatMult_SeqAIJ_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 7029123193aSHong Zhang { 7039123193aSHong Zhang PetscErrorCode ierr; 7049123193aSHong Zhang 7059123193aSHong Zhang PetscFunctionBegin; 7069123193aSHong Zhang if (scall == MAT_INITIAL_MATRIX){ 7079123193aSHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqDense(A,B,fill,C);CHKERRQ(ierr); 7089123193aSHong Zhang } 7099123193aSHong Zhang ierr = MatMatMultNumeric_SeqAIJ_SeqDense(A,B,*C);CHKERRQ(ierr); 7109123193aSHong Zhang PetscFunctionReturn(0); 7119123193aSHong Zhang } 712ec01deb9SMatthew Knepley EXTERN_C_END 713edd81eecSMatthew Knepley 7149123193aSHong Zhang #undef __FUNCT__ 7159123193aSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqDense" 7169123193aSHong Zhang PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C) 7179123193aSHong Zhang { 7189123193aSHong Zhang PetscErrorCode ierr; 7199123193aSHong Zhang 7209123193aSHong Zhang PetscFunctionBegin; 7215a586d82SBarry Smith ierr = MatMatMultSymbolic_SeqDense_SeqDense(A,B,0.0,C);CHKERRQ(ierr); 7228cdbd757SHong Zhang (*C)->ops->matmult = MatMatMult_SeqAIJ_SeqDense; 7239123193aSHong Zhang PetscFunctionReturn(0); 7249123193aSHong Zhang } 7259123193aSHong Zhang 7269123193aSHong Zhang #undef __FUNCT__ 7279123193aSHong Zhang #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqDense" 7289123193aSHong Zhang PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqDense(Mat A,Mat B,Mat C) 7299123193aSHong Zhang { 730f32d5d43SBarry Smith Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data; 7319123193aSHong Zhang PetscErrorCode ierr; 732dd6ea824SBarry Smith PetscScalar *b,*c,r1,r2,r3,r4,*b1,*b2,*b3,*b4; 733dd6ea824SBarry Smith MatScalar *aa; 734d0f46423SBarry Smith PetscInt cm=C->rmap->n, cn=B->cmap->n, bm=B->rmap->n, col, i,j,n,*aj, am = A->rmap->n; 735f32d5d43SBarry Smith PetscInt am2 = 2*am, am3 = 3*am, bm4 = 4*bm,colam; 7369123193aSHong Zhang 7379123193aSHong Zhang PetscFunctionBegin; 738f32d5d43SBarry Smith if (!cm || !cn) PetscFunctionReturn(0); 739e32f2f54SBarry Smith 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); 740e32f2f54SBarry Smith if (A->rmap->n != C->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Number rows in C %D not equal rows in A %D\n",C->rmap->n,A->rmap->n); 741e32f2f54SBarry Smith if (B->cmap->n != C->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Number columns in B %D not equal columns in C %D\n",B->cmap->n,C->cmap->n); 742f32d5d43SBarry Smith ierr = MatGetArray(B,&b);CHKERRQ(ierr); 743f32d5d43SBarry Smith ierr = MatGetArray(C,&c);CHKERRQ(ierr); 744f32d5d43SBarry Smith b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm; 745f32d5d43SBarry Smith for (col=0; col<cn-4; col += 4){ /* over columns of C */ 746f32d5d43SBarry Smith colam = col*am; 747f32d5d43SBarry Smith for (i=0; i<am; i++) { /* over rows of C in those columns */ 748f32d5d43SBarry Smith r1 = r2 = r3 = r4 = 0.0; 749f32d5d43SBarry Smith n = a->i[i+1] - a->i[i]; 750f32d5d43SBarry Smith aj = a->j + a->i[i]; 751f32d5d43SBarry Smith aa = a->a + a->i[i]; 752f32d5d43SBarry Smith for (j=0; j<n; j++) { 753f32d5d43SBarry Smith r1 += (*aa)*b1[*aj]; 754f32d5d43SBarry Smith r2 += (*aa)*b2[*aj]; 755f32d5d43SBarry Smith r3 += (*aa)*b3[*aj]; 756f32d5d43SBarry Smith r4 += (*aa++)*b4[*aj++]; 7579123193aSHong Zhang } 758f32d5d43SBarry Smith c[colam + i] = r1; 759f32d5d43SBarry Smith c[colam + am + i] = r2; 760f32d5d43SBarry Smith c[colam + am2 + i] = r3; 761f32d5d43SBarry Smith c[colam + am3 + i] = r4; 762f32d5d43SBarry Smith } 763f32d5d43SBarry Smith b1 += bm4; 764f32d5d43SBarry Smith b2 += bm4; 765f32d5d43SBarry Smith b3 += bm4; 766f32d5d43SBarry Smith b4 += bm4; 767f32d5d43SBarry Smith } 768f32d5d43SBarry Smith for (;col<cn; col++){ /* over extra columns of C */ 769f32d5d43SBarry Smith for (i=0; i<am; i++) { /* over rows of C in those columns */ 770f32d5d43SBarry Smith r1 = 0.0; 771f32d5d43SBarry Smith n = a->i[i+1] - a->i[i]; 772f32d5d43SBarry Smith aj = a->j + a->i[i]; 773f32d5d43SBarry Smith aa = a->a + a->i[i]; 774f32d5d43SBarry Smith 775f32d5d43SBarry Smith for (j=0; j<n; j++) { 776f32d5d43SBarry Smith r1 += (*aa++)*b1[*aj++]; 777f32d5d43SBarry Smith } 778f32d5d43SBarry Smith c[col*am + i] = r1; 779f32d5d43SBarry Smith } 780f32d5d43SBarry Smith b1 += bm; 781f32d5d43SBarry Smith } 782dc0b31edSSatish Balay ierr = PetscLogFlops(cn*(2.0*a->nz));CHKERRQ(ierr); 783f32d5d43SBarry Smith ierr = MatRestoreArray(B,&b);CHKERRQ(ierr); 784f32d5d43SBarry Smith ierr = MatRestoreArray(C,&c);CHKERRQ(ierr); 785f32d5d43SBarry Smith ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 786f32d5d43SBarry Smith ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 787f32d5d43SBarry Smith PetscFunctionReturn(0); 788f32d5d43SBarry Smith } 789f32d5d43SBarry Smith 790f32d5d43SBarry Smith /* 7914324174eSBarry Smith Note very similar to MatMult_SeqAIJ(), should generate both codes from same base 792f32d5d43SBarry Smith */ 793f32d5d43SBarry Smith #undef __FUNCT__ 794f32d5d43SBarry Smith #define __FUNCT__ "MatMatMultNumericAdd_SeqAIJ_SeqDense" 795f32d5d43SBarry Smith PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat A,Mat B,Mat C) 796f32d5d43SBarry Smith { 797f32d5d43SBarry Smith Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data; 798f32d5d43SBarry Smith PetscErrorCode ierr; 799dd6ea824SBarry Smith PetscScalar *b,*c,r1,r2,r3,r4,*b1,*b2,*b3,*b4; 800dd6ea824SBarry Smith MatScalar *aa; 801d0f46423SBarry Smith PetscInt cm=C->rmap->n, cn=B->cmap->n, bm=B->rmap->n, col, i,j,n,*aj, am = A->rmap->n,*ii,arm; 8024324174eSBarry Smith PetscInt am2 = 2*am, am3 = 3*am, bm4 = 4*bm,colam,*ridx; 803f32d5d43SBarry Smith 804f32d5d43SBarry Smith PetscFunctionBegin; 805f32d5d43SBarry Smith if (!cm || !cn) PetscFunctionReturn(0); 806f32d5d43SBarry Smith ierr = MatGetArray(B,&b);CHKERRQ(ierr); 807f32d5d43SBarry Smith ierr = MatGetArray(C,&c);CHKERRQ(ierr); 808f32d5d43SBarry Smith b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm; 8094324174eSBarry Smith 8104324174eSBarry Smith if (a->compressedrow.use){ /* use compressed row format */ 8114324174eSBarry Smith for (col=0; col<cn-4; col += 4){ /* over columns of C */ 8124324174eSBarry Smith colam = col*am; 8134324174eSBarry Smith arm = a->compressedrow.nrows; 8144324174eSBarry Smith ii = a->compressedrow.i; 8154324174eSBarry Smith ridx = a->compressedrow.rindex; 8164324174eSBarry Smith for (i=0; i<arm; i++) { /* over rows of C in those columns */ 8174324174eSBarry Smith r1 = r2 = r3 = r4 = 0.0; 8184324174eSBarry Smith n = ii[i+1] - ii[i]; 8194324174eSBarry Smith aj = a->j + ii[i]; 8204324174eSBarry Smith aa = a->a + ii[i]; 8214324174eSBarry Smith for (j=0; j<n; j++) { 8224324174eSBarry Smith r1 += (*aa)*b1[*aj]; 8234324174eSBarry Smith r2 += (*aa)*b2[*aj]; 8244324174eSBarry Smith r3 += (*aa)*b3[*aj]; 8254324174eSBarry Smith r4 += (*aa++)*b4[*aj++]; 8264324174eSBarry Smith } 8274324174eSBarry Smith c[colam + ridx[i]] += r1; 8284324174eSBarry Smith c[colam + am + ridx[i]] += r2; 8294324174eSBarry Smith c[colam + am2 + ridx[i]] += r3; 8304324174eSBarry Smith c[colam + am3 + ridx[i]] += r4; 8314324174eSBarry Smith } 8324324174eSBarry Smith b1 += bm4; 8334324174eSBarry Smith b2 += bm4; 8344324174eSBarry Smith b3 += bm4; 8354324174eSBarry Smith b4 += bm4; 8364324174eSBarry Smith } 8374324174eSBarry Smith for (;col<cn; col++){ /* over extra columns of C */ 8384324174eSBarry Smith colam = col*am; 8394324174eSBarry Smith arm = a->compressedrow.nrows; 8404324174eSBarry Smith ii = a->compressedrow.i; 8414324174eSBarry Smith ridx = a->compressedrow.rindex; 8424324174eSBarry Smith for (i=0; i<arm; i++) { /* over rows of C in those columns */ 8434324174eSBarry Smith r1 = 0.0; 8444324174eSBarry Smith n = ii[i+1] - ii[i]; 8454324174eSBarry Smith aj = a->j + ii[i]; 8464324174eSBarry Smith aa = a->a + ii[i]; 8474324174eSBarry Smith 8484324174eSBarry Smith for (j=0; j<n; j++) { 8494324174eSBarry Smith r1 += (*aa++)*b1[*aj++]; 8504324174eSBarry Smith } 8514324174eSBarry Smith c[col*am + ridx[i]] += r1; 8524324174eSBarry Smith } 8534324174eSBarry Smith b1 += bm; 8544324174eSBarry Smith } 8554324174eSBarry Smith } else { 856f32d5d43SBarry Smith for (col=0; col<cn-4; col += 4){ /* over columns of C */ 857f32d5d43SBarry Smith colam = col*am; 858f32d5d43SBarry Smith for (i=0; i<am; i++) { /* over rows of C in those columns */ 859f32d5d43SBarry Smith r1 = r2 = r3 = r4 = 0.0; 860f32d5d43SBarry Smith n = a->i[i+1] - a->i[i]; 861f32d5d43SBarry Smith aj = a->j + a->i[i]; 862f32d5d43SBarry Smith aa = a->a + a->i[i]; 863f32d5d43SBarry Smith for (j=0; j<n; j++) { 864f32d5d43SBarry Smith r1 += (*aa)*b1[*aj]; 865f32d5d43SBarry Smith r2 += (*aa)*b2[*aj]; 866f32d5d43SBarry Smith r3 += (*aa)*b3[*aj]; 867f32d5d43SBarry Smith r4 += (*aa++)*b4[*aj++]; 868f32d5d43SBarry Smith } 869f32d5d43SBarry Smith c[colam + i] += r1; 870f32d5d43SBarry Smith c[colam + am + i] += r2; 871f32d5d43SBarry Smith c[colam + am2 + i] += r3; 872f32d5d43SBarry Smith c[colam + am3 + i] += r4; 873f32d5d43SBarry Smith } 874f32d5d43SBarry Smith b1 += bm4; 875f32d5d43SBarry Smith b2 += bm4; 876f32d5d43SBarry Smith b3 += bm4; 877f32d5d43SBarry Smith b4 += bm4; 878f32d5d43SBarry Smith } 879f32d5d43SBarry Smith for (;col<cn; col++){ /* over extra columns of C */ 880f32d5d43SBarry Smith for (i=0; i<am; i++) { /* over rows of C in those columns */ 881f32d5d43SBarry Smith r1 = 0.0; 882f32d5d43SBarry Smith n = a->i[i+1] - a->i[i]; 883f32d5d43SBarry Smith aj = a->j + a->i[i]; 884f32d5d43SBarry Smith aa = a->a + a->i[i]; 885f32d5d43SBarry Smith 886f32d5d43SBarry Smith for (j=0; j<n; j++) { 887f32d5d43SBarry Smith r1 += (*aa++)*b1[*aj++]; 888f32d5d43SBarry Smith } 889f32d5d43SBarry Smith c[col*am + i] += r1; 890f32d5d43SBarry Smith } 891f32d5d43SBarry Smith b1 += bm; 892f32d5d43SBarry Smith } 8934324174eSBarry Smith } 894dc0b31edSSatish Balay ierr = PetscLogFlops(cn*2.0*a->nz);CHKERRQ(ierr); 895f32d5d43SBarry Smith ierr = MatRestoreArray(B,&b);CHKERRQ(ierr); 896f32d5d43SBarry Smith ierr = MatRestoreArray(C,&c);CHKERRQ(ierr); 8979123193aSHong Zhang PetscFunctionReturn(0); 8989123193aSHong Zhang } 899b1683b59SHong Zhang 900b1683b59SHong Zhang #undef __FUNCT__ 901b9af6bddSHong Zhang #define __FUNCT__ "MatTransColoringApplySpToDen_SeqAIJ" 902b9af6bddSHong Zhang PetscErrorCode MatTransColoringApplySpToDen_SeqAIJ(MatTransposeColoring coloring,Mat B,Mat Btdense) 903c8db22f6SHong Zhang { 904c8db22f6SHong Zhang PetscErrorCode ierr; 9052f5f1f90SHong Zhang Mat_SeqAIJ *b = (Mat_SeqAIJ*)B->data; 9062f5f1f90SHong Zhang Mat_SeqDense *btdense = (Mat_SeqDense*)Btdense->data; 9072f5f1f90SHong Zhang PetscInt *bi=b->i,*bj=b->j; 9082f5f1f90SHong Zhang PetscInt m=Btdense->rmap->n,n=Btdense->cmap->n,j,k,l,col,anz,*btcol,brow,ncolumns; 9092f5f1f90SHong Zhang MatScalar *btval,*btval_den,*ba=b->a; 9102f5f1f90SHong Zhang PetscInt *columns=coloring->columns,*colorforcol=coloring->colorforcol,ncolors=coloring->ncolors; 911c8db22f6SHong Zhang 912c8db22f6SHong Zhang PetscFunctionBegin; 9132f5f1f90SHong Zhang btval_den=btdense->v; 9142f5f1f90SHong Zhang ierr = PetscMemzero(btval_den,(m*n)*sizeof(MatScalar));CHKERRQ(ierr); 9152f5f1f90SHong Zhang for (k=0; k<ncolors; k++) { 9162f5f1f90SHong Zhang ncolumns = coloring->ncolumns[k]; 9172f5f1f90SHong Zhang for (l=0; l<ncolumns; l++) { /* insert a row of B to a column of Btdense */ 918d2853540SHong Zhang col = *(columns + colorforcol[k] + l); 9192f5f1f90SHong Zhang btcol = bj + bi[col]; 9202f5f1f90SHong Zhang btval = ba + bi[col]; 9212f5f1f90SHong Zhang anz = bi[col+1] - bi[col]; 922d2853540SHong Zhang for (j=0; j<anz; j++){ 9232f5f1f90SHong Zhang brow = btcol[j]; 9242f5f1f90SHong Zhang btval_den[brow] = btval[j]; 925c8db22f6SHong Zhang } 926c8db22f6SHong Zhang } 9272f5f1f90SHong Zhang btval_den += m; 928c8db22f6SHong Zhang } 929c8db22f6SHong Zhang PetscFunctionReturn(0); 930c8db22f6SHong Zhang } 931c8db22f6SHong Zhang 932c8db22f6SHong Zhang #undef __FUNCT__ 933b9af6bddSHong Zhang #define __FUNCT__ "MatTransColoringApplyDenToSp_SeqAIJ" 934b9af6bddSHong Zhang PetscErrorCode MatTransColoringApplyDenToSp_SeqAIJ(MatTransposeColoring matcoloring,Mat Cden,Mat Csp) 9358972f759SHong Zhang { 9368972f759SHong Zhang PetscErrorCode ierr; 937b2d2b04fSHong Zhang Mat_SeqAIJ *csp = (Mat_SeqAIJ*)Csp->data; 9382f5f1f90SHong Zhang PetscInt k,l,*row,*idx,m,ncolors=matcoloring->ncolors,nrows; 939b2d2b04fSHong Zhang PetscScalar *ca_den,*cp_den,*ca=csp->a; 9402f5f1f90SHong Zhang PetscInt *rows=matcoloring->rows,*spidx=matcoloring->columnsforspidx,*colorforrow=matcoloring->colorforrow; 9418972f759SHong Zhang 9428972f759SHong Zhang PetscFunctionBegin; 9438972f759SHong Zhang ierr = MatGetLocalSize(Csp,&m,PETSC_NULL);CHKERRQ(ierr); 944b2d2b04fSHong Zhang ierr = MatGetArray(Cden,&ca_den);CHKERRQ(ierr); 945b2d2b04fSHong Zhang cp_den = ca_den; 9462f5f1f90SHong Zhang for (k=0; k<ncolors; k++) { 9472f5f1f90SHong Zhang nrows = matcoloring->nrows[k]; 9482f5f1f90SHong Zhang row = rows + colorforrow[k]; 9492f5f1f90SHong Zhang idx = spidx + colorforrow[k]; 9502f5f1f90SHong Zhang for (l=0; l<nrows; l++){ 9512f5f1f90SHong Zhang ca[idx[l]] = cp_den[row[l]]; 952b2d2b04fSHong Zhang } 953b2d2b04fSHong Zhang cp_den += m; 954b2d2b04fSHong Zhang } 955b2d2b04fSHong Zhang ierr = MatRestoreArray(Cden,&ca_den);CHKERRQ(ierr); 9568972f759SHong Zhang PetscFunctionReturn(0); 9578972f759SHong Zhang } 9588972f759SHong Zhang 959*e99be685SHong Zhang /* 960*e99be685SHong Zhang MatGetColumnIJ_SeqAIJ_Color() and MatRestoreColumnIJ_SeqAIJ_Color() are customized from 961*e99be685SHong Zhang MatGetColumnIJ_SeqAIJ() and MatRestoreColumnIJ_SeqAIJ() by adding an output 962*e99be685SHong Zhang spidx[], index of a->j, to be used for setting 'columnsforspidx' in MatTransposeColoringCreate_SeqAIJ(). 963*e99be685SHong Zhang */ 964*e99be685SHong Zhang #undef __FUNCT__ 965*e99be685SHong Zhang #define __FUNCT__ "MatGetColumnIJ_SeqAIJ_Color" 966*e99be685SHong Zhang PetscErrorCode MatGetColumnIJ_SeqAIJ_Color(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool inodecompressed,PetscInt *nn,PetscInt *ia[],PetscInt *ja[],PetscInt *spidx[],PetscBool *done) 967*e99be685SHong Zhang { 968*e99be685SHong Zhang Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 969*e99be685SHong Zhang PetscErrorCode ierr; 970*e99be685SHong Zhang PetscInt i,*collengths,*cia,*cja,n = A->cmap->n,m = A->rmap->n; 971*e99be685SHong Zhang PetscInt nz = a->i[m],row,*jj,mr,col; 972*e99be685SHong Zhang PetscInt *cspidx; 973*e99be685SHong Zhang 974*e99be685SHong Zhang PetscFunctionBegin; 975*e99be685SHong Zhang *nn = n; 976*e99be685SHong Zhang if (!ia) PetscFunctionReturn(0); 977*e99be685SHong Zhang if (symmetric) { 978*e99be685SHong Zhang SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"MatGetColumnIJ_SeqAIJ_Color() not supported for the case symmetric"); 979*e99be685SHong Zhang ierr = MatToSymmetricIJ_SeqAIJ(A->rmap->n,a->i,a->j,0,oshift,ia,ja);CHKERRQ(ierr); 980*e99be685SHong Zhang } else { 981*e99be685SHong Zhang ierr = PetscMalloc((n+1)*sizeof(PetscInt),&collengths);CHKERRQ(ierr); 982*e99be685SHong Zhang ierr = PetscMemzero(collengths,n*sizeof(PetscInt));CHKERRQ(ierr); 983*e99be685SHong Zhang ierr = PetscMalloc((n+1)*sizeof(PetscInt),&cia);CHKERRQ(ierr); 984*e99be685SHong Zhang ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&cja);CHKERRQ(ierr); 985*e99be685SHong Zhang ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&cspidx);CHKERRQ(ierr); 986*e99be685SHong Zhang jj = a->j; 987*e99be685SHong Zhang for (i=0; i<nz; i++) { 988*e99be685SHong Zhang collengths[jj[i]]++; 989*e99be685SHong Zhang } 990*e99be685SHong Zhang cia[0] = oshift; 991*e99be685SHong Zhang for (i=0; i<n; i++) { 992*e99be685SHong Zhang cia[i+1] = cia[i] + collengths[i]; 993*e99be685SHong Zhang } 994*e99be685SHong Zhang ierr = PetscMemzero(collengths,n*sizeof(PetscInt));CHKERRQ(ierr); 995*e99be685SHong Zhang jj = a->j; 996*e99be685SHong Zhang for (row=0; row<m; row++) { 997*e99be685SHong Zhang mr = a->i[row+1] - a->i[row]; 998*e99be685SHong Zhang for (i=0; i<mr; i++) { 999*e99be685SHong Zhang col = *jj++; 1000*e99be685SHong Zhang cspidx[cia[col] + collengths[col] - oshift] = a->i[row] + i; /* index of a->j */ 1001*e99be685SHong Zhang cja[cia[col] + collengths[col]++ - oshift] = row + oshift; 1002*e99be685SHong Zhang } 1003*e99be685SHong Zhang } 1004*e99be685SHong Zhang ierr = PetscFree(collengths);CHKERRQ(ierr); 1005*e99be685SHong Zhang *ia = cia; *ja = cja; 1006*e99be685SHong Zhang *spidx = cspidx; 1007*e99be685SHong Zhang } 1008*e99be685SHong Zhang PetscFunctionReturn(0); 1009*e99be685SHong Zhang } 1010*e99be685SHong Zhang 1011*e99be685SHong Zhang #undef __FUNCT__ 1012*e99be685SHong Zhang #define __FUNCT__ "MatRestoreColumnIJ_SeqAIJ_Color" 1013*e99be685SHong Zhang PetscErrorCode MatRestoreColumnIJ_SeqAIJ_Color(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool inodecompressed,PetscInt *n,PetscInt *ia[],PetscInt *ja[],PetscInt *spidx[],PetscBool *done) 1014*e99be685SHong Zhang { 1015*e99be685SHong Zhang PetscErrorCode ierr; 1016*e99be685SHong Zhang 1017*e99be685SHong Zhang PetscFunctionBegin; 1018*e99be685SHong Zhang if (!ia) PetscFunctionReturn(0); 1019*e99be685SHong Zhang 1020*e99be685SHong Zhang ierr = PetscFree(*ia);CHKERRQ(ierr); 1021*e99be685SHong Zhang ierr = PetscFree(*ja);CHKERRQ(ierr); 1022*e99be685SHong Zhang ierr = PetscFree(*spidx);CHKERRQ(ierr); 1023*e99be685SHong Zhang PetscFunctionReturn(0); 1024*e99be685SHong Zhang } 1025*e99be685SHong Zhang 10268972f759SHong Zhang #undef __FUNCT__ 1027b9af6bddSHong Zhang #define __FUNCT__ "MatTransposeColoringCreate_SeqAIJ" 1028b9af6bddSHong Zhang PetscErrorCode MatTransposeColoringCreate_SeqAIJ(Mat mat,ISColoring iscoloring,MatTransposeColoring c) 1029b1683b59SHong Zhang { 1030b1683b59SHong Zhang PetscErrorCode ierr; 1031d2853540SHong Zhang PetscInt i,n,nrows,N,j,k,m,*row_idx,*ci,*cj,ncols,col,cm; 1032b1683b59SHong Zhang const PetscInt *is; 1033b28d80bdSHong Zhang PetscInt nis = iscoloring->n,*rowhit,bs = 1; 1034b1683b59SHong Zhang IS *isa; 1035b1683b59SHong Zhang PetscBool done; 1036b1683b59SHong Zhang PetscBool flg1,flg2; 1037b28d80bdSHong Zhang Mat_SeqAIJ *csp = (Mat_SeqAIJ*)mat->data; 1038*e99be685SHong Zhang PetscInt *colorforrow,*rows,*rows_i,*columnsforspidx,*columnsforspidx_i,*idxhit,*spidx; 1039d2853540SHong Zhang PetscInt *colorforcol,*columns,*columns_i; 1040b1683b59SHong Zhang 1041b1683b59SHong Zhang PetscFunctionBegin; 1042b1683b59SHong Zhang ierr = ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);CHKERRQ(ierr); 1043*e99be685SHong Zhang 1044b1683b59SHong Zhang /* this is ugly way to get blocksize but cannot call MatGetBlockSize() because AIJ can have bs > 1 */ 1045b1683b59SHong Zhang ierr = PetscTypeCompare((PetscObject)mat,MATSEQBAIJ,&flg1);CHKERRQ(ierr); 1046b1683b59SHong Zhang ierr = PetscTypeCompare((PetscObject)mat,MATMPIBAIJ,&flg2);CHKERRQ(ierr); 1047b1683b59SHong Zhang if (flg1 || flg2) { 1048b1683b59SHong Zhang ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr); 1049b1683b59SHong Zhang } 1050b1683b59SHong Zhang 1051b1683b59SHong Zhang N = mat->cmap->N/bs; 1052b1683b59SHong Zhang c->M = mat->rmap->N/bs; /* set total rows, columns and local rows */ 1053b1683b59SHong Zhang c->N = mat->cmap->N/bs; 1054b1683b59SHong Zhang c->m = mat->rmap->N/bs; 1055b1683b59SHong Zhang c->rstart = 0; 1056b1683b59SHong Zhang 1057b1683b59SHong Zhang c->ncolors = nis; 1058b1683b59SHong Zhang ierr = PetscMalloc(nis*sizeof(PetscInt),&c->ncolumns);CHKERRQ(ierr); 1059b28d80bdSHong Zhang ierr = PetscMalloc(nis*sizeof(PetscInt),&c->nrows);CHKERRQ(ierr); 1060d2853540SHong Zhang ierr = PetscMalloc2(csp->nz+1,PetscInt,&rows,csp->nz+1,PetscInt,&columnsforspidx);CHKERRQ(ierr); 1061d2853540SHong Zhang ierr = PetscMalloc((nis+1)*sizeof(PetscInt),&colorforrow);CHKERRQ(ierr); 1062d2853540SHong Zhang colorforrow[0] = 0; 1063d2853540SHong Zhang rows_i = rows; 1064d2853540SHong Zhang columnsforspidx_i = columnsforspidx; 1065b1683b59SHong Zhang 1066d2853540SHong Zhang ierr = PetscMalloc((nis+1)*sizeof(PetscInt),&colorforcol);CHKERRQ(ierr); 1067d2853540SHong Zhang ierr = PetscMalloc((N+1)*sizeof(PetscInt),&columns);CHKERRQ(ierr); 1068d2853540SHong Zhang colorforcol[0] = 0; 1069d2853540SHong Zhang columns_i = columns; 1070d2853540SHong Zhang 1071*e99be685SHong Zhang ierr = MatGetColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,&done);CHKERRQ(ierr); /* column-wise storage of mat */ 1072b1683b59SHong Zhang if (!done) SETERRQ1(((PetscObject)mat)->comm,PETSC_ERR_SUP,"MatGetColumnIJ() not supported for matrix type %s",((PetscObject)mat)->type_name); 1073b1683b59SHong Zhang 1074b28d80bdSHong Zhang cm = c->m; 1075b28d80bdSHong Zhang ierr = PetscMalloc((cm+1)*sizeof(PetscInt),&rowhit);CHKERRQ(ierr); 1076*e99be685SHong Zhang ierr = PetscMalloc((cm+1)*sizeof(PetscInt),&idxhit);CHKERRQ(ierr); 1077b1683b59SHong Zhang for (i=0; i<nis; i++) { 1078b1683b59SHong Zhang ierr = ISGetLocalSize(isa[i],&n);CHKERRQ(ierr); 1079b1683b59SHong Zhang ierr = ISGetIndices(isa[i],&is);CHKERRQ(ierr); 1080b1683b59SHong Zhang c->ncolumns[i] = n; 1081b1683b59SHong Zhang if (n) { 1082d2853540SHong Zhang ierr = PetscMemcpy(columns_i,is,n*sizeof(PetscInt));CHKERRQ(ierr); 1083b1683b59SHong Zhang } 1084d2853540SHong Zhang colorforcol[i+1] = colorforcol[i] + n; 1085d2853540SHong Zhang columns_i += n; 1086b1683b59SHong Zhang 1087b1683b59SHong Zhang /* fast, crude version requires O(N*N) work */ 1088e8354b3eSHong Zhang ierr = PetscMemzero(rowhit,cm*sizeof(PetscInt));CHKERRQ(ierr); 1089*e99be685SHong Zhang 1090b1683b59SHong Zhang /* loop over columns*/ 1091b1683b59SHong Zhang for (j=0; j<n; j++) { 1092b1683b59SHong Zhang col = is[j]; 1093d2853540SHong Zhang row_idx = cj + ci[col]; 1094b1683b59SHong Zhang m = ci[col+1] - ci[col]; 1095b1683b59SHong Zhang /* loop over columns marking them in rowhit */ 1096b1683b59SHong Zhang for (k=0; k<m; k++) { 1097*e99be685SHong Zhang idxhit[*row_idx] = spidx[ci[col] + k]; 1098d2853540SHong Zhang rowhit[*row_idx++] = col + 1; 1099b1683b59SHong Zhang } 1100b1683b59SHong Zhang } 1101b1683b59SHong Zhang /* count the number of hits */ 1102b1683b59SHong Zhang nrows = 0; 1103e8354b3eSHong Zhang for (j=0; j<cm; j++) { 1104b1683b59SHong Zhang if (rowhit[j]) nrows++; 1105b1683b59SHong Zhang } 1106b1683b59SHong Zhang c->nrows[i] = nrows; 1107d2853540SHong Zhang colorforrow[i+1] = colorforrow[i] + nrows; 1108d2853540SHong Zhang 1109b1683b59SHong Zhang nrows = 0; 1110e8354b3eSHong Zhang for (j=0; j<cm; j++) { 1111b1683b59SHong Zhang if (rowhit[j]) { 1112d2853540SHong Zhang rows_i[nrows] = j; 1113*e99be685SHong Zhang columnsforspidx_i[nrows] = idxhit[j]; 1114b1683b59SHong Zhang nrows++; 1115b1683b59SHong Zhang } 1116b1683b59SHong Zhang } 1117b1683b59SHong Zhang ierr = ISRestoreIndices(isa[i],&is);CHKERRQ(ierr); 1118d2853540SHong Zhang rows_i += nrows; columnsforspidx_i += nrows; 1119b1683b59SHong Zhang } 1120*e99be685SHong Zhang ierr = MatRestoreColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,&done);CHKERRQ(ierr); 1121b28d80bdSHong Zhang ierr = PetscFree(rowhit);CHKERRQ(ierr); 1122b1683b59SHong Zhang ierr = ISColoringRestoreIS(iscoloring,&isa);CHKERRQ(ierr); 1123d2853540SHong Zhang #if defined(PETSC_USE_DEBUG) 1124d2853540SHong Zhang if (csp->nz != colorforrow[nis]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"csp->nz %d != colorforrow[nis] %d",csp->nz,colorforrow[nis]); 1125d2853540SHong Zhang #endif 1126b28d80bdSHong Zhang 1127d2853540SHong Zhang c->colorforrow = colorforrow; 1128d2853540SHong Zhang c->rows = rows; 1129d2853540SHong Zhang c->columnsforspidx = columnsforspidx; 1130d2853540SHong Zhang c->colorforcol = colorforcol; 1131d2853540SHong Zhang c->columns = columns; 1132f94ccd6cSHong Zhang ierr = PetscFree(idxhit);CHKERRQ(ierr); 1133b1683b59SHong Zhang PetscFunctionReturn(0); 1134b1683b59SHong Zhang } 1135