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){ 2126be0446SHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr); 2226be0446SHong Zhang } 2326be0446SHong Zhang ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr); 245c66b693SKris Buschelman PetscFunctionReturn(0); 255c66b693SKris Buschelman } 26ec01deb9SMatthew Knepley EXTERN_C_END 271c24bd37SHong Zhang 2826be0446SHong Zhang #undef __FUNCT__ 2926be0446SHong Zhang #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ" 30dfbe8321SBarry Smith PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C) 3158c24d83SHong Zhang { 32dfbe8321SBarry Smith PetscErrorCode ierr; 33a1a86e44SBarry Smith PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 3458c24d83SHong Zhang Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c; 3538baddfdSBarry Smith PetscInt *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci,*cj; 36d0f46423SBarry Smith PetscInt am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N; 3738baddfdSBarry Smith PetscInt i,j,anzi,brow,bnzj,cnzi,nlnk,*lnk,nspacedouble=0; 3858c24d83SHong Zhang MatScalar *ca; 39be0fcf8dSHong Zhang PetscBT lnkbt; 407212da7cSHong Zhang PetscReal afill; 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; 612d09714cSHong Zhang aj = a->j + ai[i]; 622d09714cSHong Zhang while (j){/* assume cols are almost in increasing order, starting from its end saves computation */ 632d09714cSHong Zhang j--; 642d09714cSHong 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); 76c5db241fSHong Zhang nspacedouble++; 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 8558c24d83SHong Zhang ci[i+1] = ci[i] + cnzi; 8658c24d83SHong Zhang } 8758c24d83SHong Zhang 8858c24d83SHong Zhang /* Column indices are in the list of free space */ 8958c24d83SHong Zhang /* Allocate space for cj, initialize cj, and */ 9058c24d83SHong Zhang /* destroy list of free space and other temporary array(s) */ 9138baddfdSBarry Smith ierr = PetscMalloc((ci[am]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr); 92a1a86e44SBarry Smith ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr); 93be0fcf8dSHong Zhang ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 9458c24d83SHong Zhang 9558c24d83SHong Zhang /* Allocate space for ca */ 9658c24d83SHong Zhang ierr = PetscMalloc((ci[am]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr); 9758c24d83SHong Zhang ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr); 9858c24d83SHong Zhang 9926be0446SHong Zhang /* put together the new symbolic matrix */ 1007adad957SLisandro Dalcin ierr = MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,am,bn,ci,cj,ca,C);CHKERRQ(ierr); 10158c24d83SHong Zhang 10258c24d83SHong Zhang /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 10358c24d83SHong Zhang /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */ 10458c24d83SHong Zhang c = (Mat_SeqAIJ *)((*C)->data); 105e6b907acSBarry Smith c->free_a = PETSC_TRUE; 106e6b907acSBarry Smith c->free_ij = PETSC_TRUE; 10758c24d83SHong Zhang c->nonew = 0; 10858c24d83SHong Zhang 1097212da7cSHong Zhang /* set MatInfo */ 1107212da7cSHong Zhang afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5; 111f2b054eeSHong Zhang if (afill < 1.0) afill = 1.0; 1127212da7cSHong Zhang c->maxnz = ci[am]; 1137212da7cSHong Zhang c->nz = ci[am]; 1147212da7cSHong Zhang (*C)->info.mallocs = nspacedouble; 1157212da7cSHong Zhang (*C)->info.fill_ratio_given = fill; 1167212da7cSHong Zhang (*C)->info.fill_ratio_needed = afill; 1177212da7cSHong Zhang 1187212da7cSHong Zhang #if defined(PETSC_USE_INFO) 1197212da7cSHong Zhang if (ci[am]) { 120f2b054eeSHong Zhang ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr); 121f2b054eeSHong Zhang ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr); 122f2b054eeSHong Zhang } else { 123f2b054eeSHong Zhang ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr); 124be0fcf8dSHong Zhang } 125f2b054eeSHong Zhang #endif 12658c24d83SHong Zhang PetscFunctionReturn(0); 12758c24d83SHong Zhang } 128d50806bdSBarry Smith 1295c66b693SKris Buschelman 13026be0446SHong Zhang #undef __FUNCT__ 13126be0446SHong Zhang #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqAIJ" 132dfbe8321SBarry Smith PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C) 133d50806bdSBarry Smith { 134dfbe8321SBarry Smith PetscErrorCode ierr; 135d13dce4bSSatish Balay PetscLogDouble flops=0.0; 136d50806bdSBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 137d50806bdSBarry Smith Mat_SeqAIJ *b = (Mat_SeqAIJ *)B->data; 138d50806bdSBarry Smith Mat_SeqAIJ *c = (Mat_SeqAIJ *)C->data; 13938baddfdSBarry Smith PetscInt *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j; 140d0f46423SBarry Smith PetscInt am=A->rmap->N,cm=C->rmap->N; 141c124e916SHong Zhang PetscInt i,j,k,anzi,bnzi,cnzi,brow,nextb; 142c124e916SHong Zhang MatScalar *aa=a->a,*ba=b->a,*baj,*ca=c->a; 143d50806bdSBarry Smith 144d50806bdSBarry Smith PetscFunctionBegin; 145c124e916SHong Zhang /* clean old values in C */ 146c124e916SHong Zhang ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr); 147d50806bdSBarry Smith /* Traverse A row-wise. */ 148d50806bdSBarry Smith /* Build the ith row in C by summing over nonzero columns in A, */ 149d50806bdSBarry Smith /* the rows of B corresponding to nonzeros of A. */ 150d50806bdSBarry Smith for (i=0;i<am;i++) { 151d50806bdSBarry Smith anzi = ai[i+1] - ai[i]; 152d50806bdSBarry Smith for (j=0;j<anzi;j++) { 153d50806bdSBarry Smith brow = *aj++; 154d50806bdSBarry Smith bnzi = bi[brow+1] - bi[brow]; 155d50806bdSBarry Smith bjj = bj + bi[brow]; 156d50806bdSBarry Smith baj = ba + bi[brow]; 157c124e916SHong Zhang nextb = 0; 158c124e916SHong Zhang for (k=0; nextb<bnzi; k++) { 159c124e916SHong Zhang if (cj[k] == bjj[nextb]){ /* ccol == bcol */ 160c124e916SHong Zhang ca[k] += (*aa)*baj[nextb++]; 161c124e916SHong Zhang } 162d50806bdSBarry Smith } 163d50806bdSBarry Smith flops += 2*bnzi; 164d50806bdSBarry Smith aa++; 165d50806bdSBarry Smith } 166d50806bdSBarry Smith cnzi = ci[i+1] - ci[i]; 167d50806bdSBarry Smith ca += cnzi; 168d50806bdSBarry Smith cj += cnzi; 169d50806bdSBarry Smith } 170716bacf3SKris Buschelman ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 171716bacf3SKris Buschelman ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 172716bacf3SKris Buschelman 173d50806bdSBarry Smith ierr = PetscLogFlops(flops);CHKERRQ(ierr); 174d50806bdSBarry Smith PetscFunctionReturn(0); 175d50806bdSBarry Smith } 176bc011b1eSHong Zhang 177bc011b1eSHong Zhang #undef __FUNCT__ 178bc011b1eSHong Zhang #define __FUNCT__ "MatMatMultTranspose_SeqAIJ_SeqAIJ" 1795df89d91SHong Zhang PetscErrorCode MatMatMultTranspose_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 1805df89d91SHong Zhang { 181bc011b1eSHong Zhang PetscErrorCode ierr; 182bc011b1eSHong Zhang 183bc011b1eSHong Zhang PetscFunctionBegin; 184bc011b1eSHong Zhang if (scall == MAT_INITIAL_MATRIX){ 185bc011b1eSHong Zhang ierr = MatMatMultTransposeSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr); 186bc011b1eSHong Zhang } 187bc011b1eSHong Zhang ierr = MatMatMultTransposeNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr); 188bc011b1eSHong Zhang PetscFunctionReturn(0); 189bc011b1eSHong Zhang } 190bc011b1eSHong Zhang 191bc011b1eSHong Zhang #undef __FUNCT__ 192bc011b1eSHong Zhang #define __FUNCT__ "MatMatMultTransposeSymbolic_SeqAIJ_SeqAIJ" 193bc011b1eSHong Zhang PetscErrorCode MatMatMultTransposeSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C) 194bc011b1eSHong Zhang { 1955df89d91SHong Zhang PetscErrorCode ierr; 19681d82fe4SHong Zhang Mat Bt; 19781d82fe4SHong Zhang PetscInt *bti,*btj; 19881d82fe4SHong Zhang 19981d82fe4SHong Zhang PetscFunctionBegin; 20081d82fe4SHong Zhang /* create symbolic Bt */ 20181d82fe4SHong Zhang ierr = MatGetSymbolicTranspose_SeqAIJ(B,&bti,&btj);CHKERRQ(ierr); 20281d82fe4SHong Zhang ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,B->cmap->n,B->rmap->n,bti,btj,PETSC_NULL,&Bt);CHKERRQ(ierr); 20381d82fe4SHong Zhang 20481d82fe4SHong Zhang /* get symbolic C=A*Bt */ 20581d82fe4SHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,Bt,fill,C);CHKERRQ(ierr); 20681d82fe4SHong Zhang 20781d82fe4SHong Zhang /* clean up */ 20881d82fe4SHong Zhang ierr = MatDestroy(&Bt);CHKERRQ(ierr); 20981d82fe4SHong Zhang ierr = MatRestoreSymbolicTranspose_SeqAIJ(B,&bti,&btj);CHKERRQ(ierr); 21081d82fe4SHong Zhang 21181d82fe4SHong Zhang #if defined(INEFFICIENT_ALGORITHM) 21281d82fe4SHong Zhang /* The algorithm below computes am*bm sparse inner-product - inefficient! It will be deleted later. */ 2131fa1209cSHong Zhang PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 2141fa1209cSHong Zhang Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c; 2151fa1209cSHong Zhang PetscInt *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*ci,*cj,*acol,*bcol; 2161fa1209cSHong Zhang PetscInt am=A->rmap->N,bm=B->rmap->N; 2171fa1209cSHong Zhang PetscInt i,j,anzi,bnzj,cnzi,nlnk,*lnk,nspacedouble=0,ka,kb,index[1]; 2181fa1209cSHong Zhang MatScalar *ca; 2191fa1209cSHong Zhang PetscBT lnkbt; 2201fa1209cSHong Zhang PetscReal afill; 2215df89d91SHong Zhang 2221fa1209cSHong Zhang /* Allocate row pointer array ci */ 2231fa1209cSHong Zhang ierr = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr); 2241fa1209cSHong Zhang ci[0] = 0; 2251fa1209cSHong Zhang 2261fa1209cSHong Zhang /* Create and initialize a linked list for C columns */ 2271fa1209cSHong Zhang nlnk = bm+1; 2281fa1209cSHong Zhang ierr = PetscLLCreate(bm,bm,nlnk,lnk,lnkbt);CHKERRQ(ierr); 2291fa1209cSHong Zhang 2301fa1209cSHong Zhang /* Initial FreeSpace with size fill*(nnz(A)+nnz(B)) */ 2311fa1209cSHong Zhang ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr); 2321fa1209cSHong Zhang current_space = free_space; 2331fa1209cSHong Zhang 2341fa1209cSHong Zhang /* Determine symbolic info for each row of the product A*B^T: */ 2351fa1209cSHong Zhang for (i=0; i<am; i++) { 2361fa1209cSHong Zhang anzi = ai[i+1] - ai[i]; 2371fa1209cSHong Zhang cnzi = 0; 2381fa1209cSHong Zhang acol = aj + ai[i]; 2391fa1209cSHong Zhang for (j=0; j<bm; j++){ 2401fa1209cSHong Zhang bnzj = bi[j+1] - bi[j]; 2411fa1209cSHong Zhang bcol= bj + bi[j]; 24281d82fe4SHong Zhang /* sparse inner-product c(i,j)=A[i,:]*B[j,:]^T */ 2431fa1209cSHong Zhang ka = 0; kb = 0; 2441fa1209cSHong Zhang while (ka < anzi && kb < bnzj){ 24581d82fe4SHong Zhang while (acol[ka] < bcol[kb] && ka < anzi) ka++; 24681d82fe4SHong Zhang if (ka == anzi) break; 24781d82fe4SHong Zhang while (acol[ka] > bcol[kb] && kb < bnzj) kb++; 24881d82fe4SHong Zhang if (kb == bnzj) break; 2491fa1209cSHong Zhang if (acol[ka] == bcol[kb]){ /* add nonzero c(i,j) to lnk */ 2501fa1209cSHong Zhang index[0] = j; 2511fa1209cSHong Zhang ierr = PetscLLAdd(1,index,bm,nlnk,lnk,lnkbt);CHKERRQ(ierr); 2521fa1209cSHong Zhang cnzi++; 2531fa1209cSHong Zhang break; 2541fa1209cSHong Zhang } 2551fa1209cSHong Zhang } 2561fa1209cSHong Zhang } 2571fa1209cSHong Zhang 2581fa1209cSHong Zhang /* If free space is not available, make more free space */ 2591fa1209cSHong Zhang /* Double the amount of total space in the list */ 2601fa1209cSHong Zhang if (current_space->local_remaining<cnzi) { 2611fa1209cSHong Zhang ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 2621fa1209cSHong Zhang nspacedouble++; 2631fa1209cSHong Zhang } 2641fa1209cSHong Zhang 2651fa1209cSHong Zhang /* Copy data into free space, then initialize lnk */ 2661fa1209cSHong Zhang ierr = PetscLLClean(bm,bm,cnzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 2671fa1209cSHong Zhang current_space->array += cnzi; 2681fa1209cSHong Zhang current_space->local_used += cnzi; 2691fa1209cSHong Zhang current_space->local_remaining -= cnzi; 2701fa1209cSHong Zhang 2711fa1209cSHong Zhang ci[i+1] = ci[i] + cnzi; 2721fa1209cSHong Zhang } 2731fa1209cSHong Zhang 2741fa1209cSHong Zhang 2751fa1209cSHong Zhang /* Column indices are in the list of free space. 2761fa1209cSHong Zhang Allocate array cj, copy column indices to cj, and destroy list of free space */ 2771fa1209cSHong Zhang ierr = PetscMalloc((ci[am]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr); 2781fa1209cSHong Zhang ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr); 2791fa1209cSHong Zhang ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 2801fa1209cSHong Zhang 2811fa1209cSHong Zhang /* Allocate space for ca */ 2821fa1209cSHong Zhang ierr = PetscMalloc((ci[am]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr); 2831fa1209cSHong Zhang ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr); 2841fa1209cSHong Zhang 2851fa1209cSHong Zhang /* put together the new symbolic matrix */ 2861fa1209cSHong Zhang ierr = MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,am,bm,ci,cj,ca,C);CHKERRQ(ierr); 2871fa1209cSHong Zhang 2881fa1209cSHong Zhang /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 2891fa1209cSHong Zhang /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */ 2901fa1209cSHong Zhang c = (Mat_SeqAIJ *)((*C)->data); 2911fa1209cSHong Zhang c->free_a = PETSC_TRUE; 2921fa1209cSHong Zhang c->free_ij = PETSC_TRUE; 2931fa1209cSHong Zhang c->nonew = 0; 2941fa1209cSHong Zhang 2951fa1209cSHong Zhang /* set MatInfo */ 2961fa1209cSHong Zhang afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5; 2971fa1209cSHong Zhang if (afill < 1.0) afill = 1.0; 2981fa1209cSHong Zhang c->maxnz = ci[am]; 2991fa1209cSHong Zhang c->nz = ci[am]; 3001fa1209cSHong Zhang (*C)->info.mallocs = nspacedouble; 3011fa1209cSHong Zhang (*C)->info.fill_ratio_given = fill; 3021fa1209cSHong Zhang (*C)->info.fill_ratio_needed = afill; 3031fa1209cSHong Zhang 3041fa1209cSHong Zhang #if defined(PETSC_USE_INFO) 3051fa1209cSHong Zhang if (ci[am]) { 3061fa1209cSHong Zhang ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr); 3071fa1209cSHong Zhang ierr = PetscInfo1((*C),"Use MatMatMultTranspose(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr); 3081fa1209cSHong Zhang } else { 3091fa1209cSHong Zhang ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr); 3101fa1209cSHong Zhang } 3111fa1209cSHong Zhang #endif 31281d82fe4SHong Zhang #endif 3135df89d91SHong Zhang PetscFunctionReturn(0); 3145df89d91SHong Zhang } 3155df89d91SHong Zhang 3166973769bSHong Zhang /* #define USE_ARRAY - for sparse dot product. Slower than !USE_ARRAY */ 3175df89d91SHong Zhang #undef __FUNCT__ 3185df89d91SHong Zhang #define __FUNCT__ "MatMatMultTransposeNumeric_SeqAIJ_SeqAIJ" 3195df89d91SHong Zhang PetscErrorCode MatMatMultTransposeNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C) 3205df89d91SHong Zhang { 3215df89d91SHong Zhang PetscErrorCode ierr; 3225df89d91SHong Zhang Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data; 323e2cac8caSJed Brown PetscInt *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,anzi,bnzj,nexta,nextb,*acol,*bcol,brow; 32481d82fe4SHong Zhang PetscInt cm=C->rmap->n,*ci=c->i,*cj=c->j,i,j,cnzi,*ccol; 3255df89d91SHong Zhang PetscLogDouble flops=0.0; 3266973769bSHong Zhang MatScalar *aa=a->a,*aval,*ba=b->a,*bval,*ca=c->a,*cval; 3276973769bSHong Zhang #if defined(USE_ARRAY) 3286973769bSHong Zhang MatScalar *spdot; 3296973769bSHong Zhang #endif 330*c8db22f6SHong Zhang PetscBool flg; 3315df89d91SHong Zhang 3325df89d91SHong Zhang PetscFunctionBegin; 333*c8db22f6SHong Zhang ierr = PetscOptionsGetBool(PETSC_NULL,"-matmulttrans_color",&flg,PETSC_NULL);CHKERRQ(ierr); 334*c8db22f6SHong Zhang 335*c8db22f6SHong Zhang if (flg){ 336*c8db22f6SHong Zhang printf("Create MatMultTransposeColoring ...\n"); 337*c8db22f6SHong Zhang /* Create MatMultTransposeColoring from symbolic C=A*B^T */ 338*c8db22f6SHong Zhang MatMultTransposeColoring matfdcoloring = 0; 339*c8db22f6SHong Zhang ISColoring iscoloring; 340*c8db22f6SHong Zhang ierr = MatGetColoring(C,MATCOLORINGSL,&iscoloring);CHKERRQ(ierr); 341*c8db22f6SHong Zhang ierr = MatMultTransposeColoringCreate(C,iscoloring,&matfdcoloring);CHKERRQ(ierr); 342*c8db22f6SHong Zhang ierr = ISColoringDestroy(&iscoloring);CHKERRQ(ierr); 343*c8db22f6SHong Zhang 344*c8db22f6SHong Zhang /* Create Bt_dense */ 345*c8db22f6SHong Zhang Mat Bt_dense; 346*c8db22f6SHong Zhang PetscInt m,n; 347*c8db22f6SHong Zhang ierr = MatCreate(PETSC_COMM_WORLD,&Bt_dense);CHKERRQ(ierr); 348*c8db22f6SHong Zhang ierr = MatSetSizes(Bt_dense,A->cmap->n,matfdcoloring->ncolors,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 349*c8db22f6SHong Zhang ierr = MatSetType(Bt_dense,MATDENSE);CHKERRQ(ierr); 350*c8db22f6SHong Zhang ierr = MatSeqDenseSetPreallocation(Bt_dense,PETSC_NULL);CHKERRQ(ierr); 351*c8db22f6SHong Zhang ierr = MatAssemblyBegin(Bt_dense,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 352*c8db22f6SHong Zhang ierr = MatAssemblyEnd(Bt_dense,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 353*c8db22f6SHong Zhang ierr = MatGetLocalSize(Bt_dense,&m,&n);CHKERRQ(ierr); 354*c8db22f6SHong Zhang printf("Bt_dense: %d,%d\n",m,n); 355*c8db22f6SHong Zhang 356*c8db22f6SHong Zhang /* Get Bt_dense by Apply MatMultTransposeColoring to B */ 357*c8db22f6SHong Zhang ierr = MatMultTransposeColoringApply(B,Bt_dense,matfdcoloring);CHKERRQ(ierr); 358*c8db22f6SHong Zhang 359*c8db22f6SHong Zhang /* C_dense = A*Bt_dense */ 360*c8db22f6SHong Zhang Mat C_dense; 361*c8db22f6SHong Zhang ierr = MatMatMult(A,Bt_dense,MAT_INITIAL_MATRIX,2.0,&C_dense); CHKERRQ(ierr); 362*c8db22f6SHong Zhang ierr = MatSetOptionsPrefix(C_dense,"C_dense_");CHKERRQ(ierr); 363*c8db22f6SHong Zhang 364*c8db22f6SHong Zhang 365*c8db22f6SHong Zhang 366*c8db22f6SHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not done yet"); 367*c8db22f6SHong Zhang PetscFunctionReturn(0); 368*c8db22f6SHong Zhang } 369*c8db22f6SHong Zhang 3706973769bSHong Zhang #if defined(USE_ARRAY) 3716973769bSHong Zhang /* allocate an array for implementing sparse inner-product */ 372e2cac8caSJed Brown ierr = PetscMalloc((A->cmap->n+1)*sizeof(MatScalar),&spdot);CHKERRQ(ierr); 373e2cac8caSJed Brown ierr = PetscMemzero(spdot,(A->cmap->n+1)*sizeof(MatScalar));CHKERRQ(ierr); 3746973769bSHong Zhang #endif 3756973769bSHong Zhang 37681d82fe4SHong Zhang /* clear old values in C */ 37781d82fe4SHong Zhang ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr); 3785df89d91SHong Zhang 3791fa1209cSHong Zhang for (i=0; i<cm; i++) { 38081d82fe4SHong Zhang anzi = ai[i+1] - ai[i]; 38181d82fe4SHong Zhang acol = aj + ai[i]; 3826973769bSHong Zhang aval = aa + ai[i]; 3831fa1209cSHong Zhang cnzi = ci[i+1] - ci[i]; 3841fa1209cSHong Zhang ccol = cj + ci[i]; 3856973769bSHong Zhang cval = ca + ci[i]; 3861fa1209cSHong Zhang for (j=0; j<cnzi; j++){ 38781d82fe4SHong Zhang brow = ccol[j]; 38881d82fe4SHong Zhang bnzj = bi[brow+1] - bi[brow]; 38981d82fe4SHong Zhang bcol = bj + bi[brow]; 3906973769bSHong Zhang bval = ba + bi[brow]; 3916973769bSHong Zhang 39281d82fe4SHong Zhang /* perform sparse inner-product c(i,j)=A[i,:]*B[j,:]^T */ 3936973769bSHong Zhang #if defined(USE_ARRAY) 3946973769bSHong Zhang /* put ba to spdot array */ 3956973769bSHong Zhang for (nextb=0; nextb<bnzj; nextb++) spdot[bcol[nextb]] = bval[nextb]; 3966973769bSHong Zhang /* c(i,j)=A[i,:]*B[j,:]^T */ 3976973769bSHong Zhang for (nexta=0; nexta<anzi; nexta++){ 3986973769bSHong Zhang cval[j] += spdot[acol[nexta]]*aval[nexta]; 3996973769bSHong Zhang } 4006973769bSHong Zhang /* zero spdot array */ 4016973769bSHong Zhang for (nextb=0; nextb<bnzj; nextb++) spdot[bcol[nextb]] = 0.0; 4026973769bSHong Zhang #else 40381d82fe4SHong Zhang nexta = 0; nextb = 0; 40481d82fe4SHong Zhang while (nexta<anzi && nextb<bnzj){ 40581d82fe4SHong Zhang while (acol[nexta] < bcol[nextb] && nexta < anzi) nexta++; 40681d82fe4SHong Zhang if (nexta == anzi) break; 40781d82fe4SHong Zhang while (acol[nexta] > bcol[nextb] && nextb < bnzj) nextb++; 40881d82fe4SHong Zhang if (nextb == bnzj) break; 40981d82fe4SHong Zhang if (acol[nexta] == bcol[nextb]){ 4106973769bSHong Zhang cval[j] += aval[nexta]*bval[nextb]; 41181d82fe4SHong Zhang nexta++; nextb++; 41281d82fe4SHong Zhang flops += 2; 4131fa1209cSHong Zhang } 4141fa1209cSHong Zhang } 4156973769bSHong Zhang #endif 41681d82fe4SHong Zhang } 41781d82fe4SHong Zhang } 41881d82fe4SHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 41981d82fe4SHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 42081d82fe4SHong Zhang ierr = PetscLogFlops(flops);CHKERRQ(ierr); 4216973769bSHong Zhang #if defined(USE_ARRAY) 4226973769bSHong Zhang ierr = PetscFree(spdot); 4236973769bSHong Zhang #endif 4245df89d91SHong Zhang PetscFunctionReturn(0); 4255df89d91SHong Zhang } 4265df89d91SHong Zhang 4275df89d91SHong Zhang #undef __FUNCT__ 4285df89d91SHong Zhang #define __FUNCT__ "MatMatTransposeMult_SeqAIJ_SeqAIJ" 4295df89d91SHong Zhang PetscErrorCode MatMatTransposeMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) { 4305df89d91SHong Zhang PetscErrorCode ierr; 4315df89d91SHong Zhang 4325df89d91SHong Zhang PetscFunctionBegin; 4335df89d91SHong Zhang if (scall == MAT_INITIAL_MATRIX){ 4345df89d91SHong Zhang ierr = MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr); 4355df89d91SHong Zhang } 4365df89d91SHong Zhang ierr = MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr); 4375df89d91SHong Zhang PetscFunctionReturn(0); 4385df89d91SHong Zhang } 4395df89d91SHong Zhang 4405df89d91SHong Zhang #undef __FUNCT__ 4415df89d91SHong Zhang #define __FUNCT__ "MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ" 4425df89d91SHong Zhang PetscErrorCode MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C) 4435df89d91SHong Zhang { 444bc011b1eSHong Zhang PetscErrorCode ierr; 445bc011b1eSHong Zhang Mat At; 44638baddfdSBarry Smith PetscInt *ati,*atj; 447bc011b1eSHong Zhang 448bc011b1eSHong Zhang PetscFunctionBegin; 449bc011b1eSHong Zhang /* create symbolic At */ 450bc011b1eSHong Zhang ierr = MatGetSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr); 451d0f46423SBarry Smith ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,A->cmap->n,A->rmap->n,ati,atj,PETSC_NULL,&At);CHKERRQ(ierr); 452bc011b1eSHong Zhang 453bc011b1eSHong Zhang /* get symbolic C=At*B */ 454bc011b1eSHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(At,B,fill,C);CHKERRQ(ierr); 455bc011b1eSHong Zhang 456bc011b1eSHong Zhang /* clean up */ 4576bf464f9SBarry Smith ierr = MatDestroy(&At);CHKERRQ(ierr); 458bc011b1eSHong Zhang ierr = MatRestoreSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr); 459bc011b1eSHong Zhang PetscFunctionReturn(0); 460bc011b1eSHong Zhang } 461bc011b1eSHong Zhang 462bc011b1eSHong Zhang #undef __FUNCT__ 4635df89d91SHong Zhang #define __FUNCT__ "MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ" 4645df89d91SHong Zhang PetscErrorCode MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C) 465bc011b1eSHong Zhang { 466bc011b1eSHong Zhang PetscErrorCode ierr; 4670fbc74f4SHong Zhang Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data; 468d0f46423SBarry Smith PetscInt am=A->rmap->n,anzi,*ai=a->i,*aj=a->j,*bi=b->i,*bj,bnzi,nextb; 469d0f46423SBarry Smith PetscInt cm=C->rmap->n,*ci=c->i,*cj=c->j,crow,*cjj,i,j,k; 470d13dce4bSSatish Balay PetscLogDouble flops=0.0; 4710fbc74f4SHong Zhang MatScalar *aa=a->a,*ba,*ca=c->a,*caj; 472bc011b1eSHong Zhang 473bc011b1eSHong Zhang PetscFunctionBegin; 474bc011b1eSHong Zhang /* clear old values in C */ 475bc011b1eSHong Zhang ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr); 476bc011b1eSHong Zhang 477bc011b1eSHong Zhang /* compute A^T*B using outer product (A^T)[:,i]*B[i,:] */ 478bc011b1eSHong Zhang for (i=0;i<am;i++) { 479bc011b1eSHong Zhang bj = b->j + bi[i]; 480bc011b1eSHong Zhang ba = b->a + bi[i]; 481bc011b1eSHong Zhang bnzi = bi[i+1] - bi[i]; 482bc011b1eSHong Zhang anzi = ai[i+1] - ai[i]; 483bc011b1eSHong Zhang for (j=0; j<anzi; j++) { 484bc011b1eSHong Zhang nextb = 0; 4850fbc74f4SHong Zhang crow = *aj++; 486bc011b1eSHong Zhang cjj = cj + ci[crow]; 487bc011b1eSHong Zhang caj = ca + ci[crow]; 488bc011b1eSHong Zhang /* perform sparse axpy operation. Note cjj includes bj. */ 489bc011b1eSHong Zhang for (k=0; nextb<bnzi; k++) { 4900fbc74f4SHong Zhang if (cjj[k] == *(bj+nextb)) { /* ccol == bcol */ 4910fbc74f4SHong Zhang caj[k] += (*aa)*(*(ba+nextb)); 492bc011b1eSHong Zhang nextb++; 493bc011b1eSHong Zhang } 494bc011b1eSHong Zhang } 495bc011b1eSHong Zhang flops += 2*bnzi; 4960fbc74f4SHong Zhang aa++; 497bc011b1eSHong Zhang } 498bc011b1eSHong Zhang } 499bc011b1eSHong Zhang 500bc011b1eSHong Zhang /* Assemble the final matrix and clean up */ 501bc011b1eSHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 502bc011b1eSHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 503bc011b1eSHong Zhang ierr = PetscLogFlops(flops);CHKERRQ(ierr); 504bc011b1eSHong Zhang PetscFunctionReturn(0); 505bc011b1eSHong Zhang } 5069123193aSHong Zhang 507ec01deb9SMatthew Knepley EXTERN_C_BEGIN 5089123193aSHong Zhang #undef __FUNCT__ 5099123193aSHong Zhang #define __FUNCT__ "MatMatMult_SeqAIJ_SeqDense" 5109123193aSHong Zhang PetscErrorCode MatMatMult_SeqAIJ_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 5119123193aSHong Zhang { 5129123193aSHong Zhang PetscErrorCode ierr; 5139123193aSHong Zhang 5149123193aSHong Zhang PetscFunctionBegin; 5159123193aSHong Zhang if (scall == MAT_INITIAL_MATRIX){ 5169123193aSHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqDense(A,B,fill,C);CHKERRQ(ierr); 5179123193aSHong Zhang } 5189123193aSHong Zhang ierr = MatMatMultNumeric_SeqAIJ_SeqDense(A,B,*C);CHKERRQ(ierr); 5199123193aSHong Zhang PetscFunctionReturn(0); 5209123193aSHong Zhang } 521ec01deb9SMatthew Knepley EXTERN_C_END 522edd81eecSMatthew Knepley 5239123193aSHong Zhang #undef __FUNCT__ 5249123193aSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqDense" 5259123193aSHong Zhang PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C) 5269123193aSHong Zhang { 5279123193aSHong Zhang PetscErrorCode ierr; 5289123193aSHong Zhang 5299123193aSHong Zhang PetscFunctionBegin; 5305a586d82SBarry Smith ierr = MatMatMultSymbolic_SeqDense_SeqDense(A,B,0.0,C);CHKERRQ(ierr); 5319123193aSHong Zhang PetscFunctionReturn(0); 5329123193aSHong Zhang } 5339123193aSHong Zhang 5349123193aSHong Zhang #undef __FUNCT__ 5359123193aSHong Zhang #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqDense" 5369123193aSHong Zhang PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqDense(Mat A,Mat B,Mat C) 5379123193aSHong Zhang { 538f32d5d43SBarry Smith Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data; 5399123193aSHong Zhang PetscErrorCode ierr; 540dd6ea824SBarry Smith PetscScalar *b,*c,r1,r2,r3,r4,*b1,*b2,*b3,*b4; 541dd6ea824SBarry Smith MatScalar *aa; 542d0f46423SBarry Smith PetscInt cm=C->rmap->n, cn=B->cmap->n, bm=B->rmap->n, col, i,j,n,*aj, am = A->rmap->n; 543f32d5d43SBarry Smith PetscInt am2 = 2*am, am3 = 3*am, bm4 = 4*bm,colam; 5449123193aSHong Zhang 5459123193aSHong Zhang PetscFunctionBegin; 546f32d5d43SBarry Smith if (!cm || !cn) PetscFunctionReturn(0); 547e32f2f54SBarry 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); 548e32f2f54SBarry 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); 549e32f2f54SBarry 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); 550f32d5d43SBarry Smith ierr = MatGetArray(B,&b);CHKERRQ(ierr); 551f32d5d43SBarry Smith ierr = MatGetArray(C,&c);CHKERRQ(ierr); 552f32d5d43SBarry Smith b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm; 553f32d5d43SBarry Smith for (col=0; col<cn-4; col += 4){ /* over columns of C */ 554f32d5d43SBarry Smith colam = col*am; 555f32d5d43SBarry Smith for (i=0; i<am; i++) { /* over rows of C in those columns */ 556f32d5d43SBarry Smith r1 = r2 = r3 = r4 = 0.0; 557f32d5d43SBarry Smith n = a->i[i+1] - a->i[i]; 558f32d5d43SBarry Smith aj = a->j + a->i[i]; 559f32d5d43SBarry Smith aa = a->a + a->i[i]; 560f32d5d43SBarry Smith for (j=0; j<n; j++) { 561f32d5d43SBarry Smith r1 += (*aa)*b1[*aj]; 562f32d5d43SBarry Smith r2 += (*aa)*b2[*aj]; 563f32d5d43SBarry Smith r3 += (*aa)*b3[*aj]; 564f32d5d43SBarry Smith r4 += (*aa++)*b4[*aj++]; 5659123193aSHong Zhang } 566f32d5d43SBarry Smith c[colam + i] = r1; 567f32d5d43SBarry Smith c[colam + am + i] = r2; 568f32d5d43SBarry Smith c[colam + am2 + i] = r3; 569f32d5d43SBarry Smith c[colam + am3 + i] = r4; 570f32d5d43SBarry Smith } 571f32d5d43SBarry Smith b1 += bm4; 572f32d5d43SBarry Smith b2 += bm4; 573f32d5d43SBarry Smith b3 += bm4; 574f32d5d43SBarry Smith b4 += bm4; 575f32d5d43SBarry Smith } 576f32d5d43SBarry Smith for (;col<cn; col++){ /* over extra columns of C */ 577f32d5d43SBarry Smith for (i=0; i<am; i++) { /* over rows of C in those columns */ 578f32d5d43SBarry Smith r1 = 0.0; 579f32d5d43SBarry Smith n = a->i[i+1] - a->i[i]; 580f32d5d43SBarry Smith aj = a->j + a->i[i]; 581f32d5d43SBarry Smith aa = a->a + a->i[i]; 582f32d5d43SBarry Smith 583f32d5d43SBarry Smith for (j=0; j<n; j++) { 584f32d5d43SBarry Smith r1 += (*aa++)*b1[*aj++]; 585f32d5d43SBarry Smith } 586f32d5d43SBarry Smith c[col*am + i] = r1; 587f32d5d43SBarry Smith } 588f32d5d43SBarry Smith b1 += bm; 589f32d5d43SBarry Smith } 590dc0b31edSSatish Balay ierr = PetscLogFlops(cn*(2.0*a->nz));CHKERRQ(ierr); 591f32d5d43SBarry Smith ierr = MatRestoreArray(B,&b);CHKERRQ(ierr); 592f32d5d43SBarry Smith ierr = MatRestoreArray(C,&c);CHKERRQ(ierr); 593f32d5d43SBarry Smith ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 594f32d5d43SBarry Smith ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 595f32d5d43SBarry Smith PetscFunctionReturn(0); 596f32d5d43SBarry Smith } 597f32d5d43SBarry Smith 598f32d5d43SBarry Smith /* 5994324174eSBarry Smith Note very similar to MatMult_SeqAIJ(), should generate both codes from same base 600f32d5d43SBarry Smith */ 601f32d5d43SBarry Smith #undef __FUNCT__ 602f32d5d43SBarry Smith #define __FUNCT__ "MatMatMultNumericAdd_SeqAIJ_SeqDense" 603f32d5d43SBarry Smith PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat A,Mat B,Mat C) 604f32d5d43SBarry Smith { 605f32d5d43SBarry Smith Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data; 606f32d5d43SBarry Smith PetscErrorCode ierr; 607dd6ea824SBarry Smith PetscScalar *b,*c,r1,r2,r3,r4,*b1,*b2,*b3,*b4; 608dd6ea824SBarry Smith MatScalar *aa; 609d0f46423SBarry 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; 6104324174eSBarry Smith PetscInt am2 = 2*am, am3 = 3*am, bm4 = 4*bm,colam,*ridx; 611f32d5d43SBarry Smith 612f32d5d43SBarry Smith PetscFunctionBegin; 613f32d5d43SBarry Smith if (!cm || !cn) PetscFunctionReturn(0); 614f32d5d43SBarry Smith ierr = MatGetArray(B,&b);CHKERRQ(ierr); 615f32d5d43SBarry Smith ierr = MatGetArray(C,&c);CHKERRQ(ierr); 616f32d5d43SBarry Smith b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm; 6174324174eSBarry Smith 6184324174eSBarry Smith if (a->compressedrow.use){ /* use compressed row format */ 6194324174eSBarry Smith for (col=0; col<cn-4; col += 4){ /* over columns of C */ 6204324174eSBarry Smith colam = col*am; 6214324174eSBarry Smith arm = a->compressedrow.nrows; 6224324174eSBarry Smith ii = a->compressedrow.i; 6234324174eSBarry Smith ridx = a->compressedrow.rindex; 6244324174eSBarry Smith for (i=0; i<arm; i++) { /* over rows of C in those columns */ 6254324174eSBarry Smith r1 = r2 = r3 = r4 = 0.0; 6264324174eSBarry Smith n = ii[i+1] - ii[i]; 6274324174eSBarry Smith aj = a->j + ii[i]; 6284324174eSBarry Smith aa = a->a + ii[i]; 6294324174eSBarry Smith for (j=0; j<n; j++) { 6304324174eSBarry Smith r1 += (*aa)*b1[*aj]; 6314324174eSBarry Smith r2 += (*aa)*b2[*aj]; 6324324174eSBarry Smith r3 += (*aa)*b3[*aj]; 6334324174eSBarry Smith r4 += (*aa++)*b4[*aj++]; 6344324174eSBarry Smith } 6354324174eSBarry Smith c[colam + ridx[i]] += r1; 6364324174eSBarry Smith c[colam + am + ridx[i]] += r2; 6374324174eSBarry Smith c[colam + am2 + ridx[i]] += r3; 6384324174eSBarry Smith c[colam + am3 + ridx[i]] += r4; 6394324174eSBarry Smith } 6404324174eSBarry Smith b1 += bm4; 6414324174eSBarry Smith b2 += bm4; 6424324174eSBarry Smith b3 += bm4; 6434324174eSBarry Smith b4 += bm4; 6444324174eSBarry Smith } 6454324174eSBarry Smith for (;col<cn; col++){ /* over extra columns of C */ 6464324174eSBarry Smith colam = col*am; 6474324174eSBarry Smith arm = a->compressedrow.nrows; 6484324174eSBarry Smith ii = a->compressedrow.i; 6494324174eSBarry Smith ridx = a->compressedrow.rindex; 6504324174eSBarry Smith for (i=0; i<arm; i++) { /* over rows of C in those columns */ 6514324174eSBarry Smith r1 = 0.0; 6524324174eSBarry Smith n = ii[i+1] - ii[i]; 6534324174eSBarry Smith aj = a->j + ii[i]; 6544324174eSBarry Smith aa = a->a + ii[i]; 6554324174eSBarry Smith 6564324174eSBarry Smith for (j=0; j<n; j++) { 6574324174eSBarry Smith r1 += (*aa++)*b1[*aj++]; 6584324174eSBarry Smith } 6594324174eSBarry Smith c[col*am + ridx[i]] += r1; 6604324174eSBarry Smith } 6614324174eSBarry Smith b1 += bm; 6624324174eSBarry Smith } 6634324174eSBarry Smith } else { 664f32d5d43SBarry Smith for (col=0; col<cn-4; col += 4){ /* over columns of C */ 665f32d5d43SBarry Smith colam = col*am; 666f32d5d43SBarry Smith for (i=0; i<am; i++) { /* over rows of C in those columns */ 667f32d5d43SBarry Smith r1 = r2 = r3 = r4 = 0.0; 668f32d5d43SBarry Smith n = a->i[i+1] - a->i[i]; 669f32d5d43SBarry Smith aj = a->j + a->i[i]; 670f32d5d43SBarry Smith aa = a->a + a->i[i]; 671f32d5d43SBarry Smith for (j=0; j<n; j++) { 672f32d5d43SBarry Smith r1 += (*aa)*b1[*aj]; 673f32d5d43SBarry Smith r2 += (*aa)*b2[*aj]; 674f32d5d43SBarry Smith r3 += (*aa)*b3[*aj]; 675f32d5d43SBarry Smith r4 += (*aa++)*b4[*aj++]; 676f32d5d43SBarry Smith } 677f32d5d43SBarry Smith c[colam + i] += r1; 678f32d5d43SBarry Smith c[colam + am + i] += r2; 679f32d5d43SBarry Smith c[colam + am2 + i] += r3; 680f32d5d43SBarry Smith c[colam + am3 + i] += r4; 681f32d5d43SBarry Smith } 682f32d5d43SBarry Smith b1 += bm4; 683f32d5d43SBarry Smith b2 += bm4; 684f32d5d43SBarry Smith b3 += bm4; 685f32d5d43SBarry Smith b4 += bm4; 686f32d5d43SBarry Smith } 687f32d5d43SBarry Smith for (;col<cn; col++){ /* over extra columns of C */ 688f32d5d43SBarry Smith for (i=0; i<am; i++) { /* over rows of C in those columns */ 689f32d5d43SBarry Smith r1 = 0.0; 690f32d5d43SBarry Smith n = a->i[i+1] - a->i[i]; 691f32d5d43SBarry Smith aj = a->j + a->i[i]; 692f32d5d43SBarry Smith aa = a->a + a->i[i]; 693f32d5d43SBarry Smith 694f32d5d43SBarry Smith for (j=0; j<n; j++) { 695f32d5d43SBarry Smith r1 += (*aa++)*b1[*aj++]; 696f32d5d43SBarry Smith } 697f32d5d43SBarry Smith c[col*am + i] += r1; 698f32d5d43SBarry Smith } 699f32d5d43SBarry Smith b1 += bm; 700f32d5d43SBarry Smith } 7014324174eSBarry Smith } 702dc0b31edSSatish Balay ierr = PetscLogFlops(cn*2.0*a->nz);CHKERRQ(ierr); 703f32d5d43SBarry Smith ierr = MatRestoreArray(B,&b);CHKERRQ(ierr); 704f32d5d43SBarry Smith ierr = MatRestoreArray(C,&c);CHKERRQ(ierr); 7059123193aSHong Zhang PetscFunctionReturn(0); 7069123193aSHong Zhang } 707b1683b59SHong Zhang 708b1683b59SHong Zhang #undef __FUNCT__ 709*c8db22f6SHong Zhang #define __FUNCT__ "MatMultTransposeColoringApply_SeqAIJ" 710*c8db22f6SHong Zhang PetscErrorCode MatMultTransposeColoringApply_SeqAIJ(Mat B,Mat Btdense,MatMultTransposeColoring coloring) 711*c8db22f6SHong Zhang { 712*c8db22f6SHong Zhang PetscErrorCode ierr; 713*c8db22f6SHong Zhang Mat_SeqAIJ *a = (Mat_SeqAIJ*)B->data; 714*c8db22f6SHong Zhang Mat_SeqDense *atdense = (Mat_SeqDense*)Btdense->data; 715*c8db22f6SHong Zhang PetscInt m=Btdense->rmap->n,n=Btdense->cmap->n,j,k,l,col,anz,*atcol,brow,bcol; 716*c8db22f6SHong Zhang MatScalar *atval,*bval; 717*c8db22f6SHong Zhang 718*c8db22f6SHong Zhang PetscFunctionBegin; 719*c8db22f6SHong Zhang ierr = PetscMemzero(atdense->v,(m*n)*sizeof(MatScalar));CHKERRQ(ierr); 720*c8db22f6SHong Zhang for (k=0; k<coloring->ncolors; k++) { 721*c8db22f6SHong Zhang for (l=0; l<coloring->ncolumns[k]; l++) { /* insert a row of B to a column of Btdense */ 722*c8db22f6SHong Zhang col = coloring->columns[k][l]; /* =row of B */ 723*c8db22f6SHong Zhang anz = a->i[col+1] - a->i[col]; 724*c8db22f6SHong Zhang for (j=0; j<anz; j++){ 725*c8db22f6SHong Zhang atcol = a->j + a->i[col]; 726*c8db22f6SHong Zhang atval = a->a + a->i[col]; 727*c8db22f6SHong Zhang bval = atdense->v; 728*c8db22f6SHong Zhang brow = atcol[j]; 729*c8db22f6SHong Zhang bcol = k; 730*c8db22f6SHong Zhang bval[bcol*m+brow] = atval[j]; 731*c8db22f6SHong Zhang } 732*c8db22f6SHong Zhang } 733*c8db22f6SHong Zhang } 734*c8db22f6SHong Zhang PetscFunctionReturn(0); 735*c8db22f6SHong Zhang } 736*c8db22f6SHong Zhang 737*c8db22f6SHong Zhang #undef __FUNCT__ 738b1683b59SHong Zhang #define __FUNCT__ "MatMultTransposeColoringCreate_SeqAIJ" 739b1683b59SHong Zhang PetscErrorCode MatMultTransposeColoringCreate_SeqAIJ(Mat mat,ISColoring iscoloring,MatMultTransposeColoring c) 740b1683b59SHong Zhang { 741b1683b59SHong Zhang PetscErrorCode ierr; 742b1683b59SHong Zhang PetscInt i,n,nrows,N,j,k,m,*rows,*ci,*cj,ncols,col; 743b1683b59SHong Zhang const PetscInt *is; 744b1683b59SHong Zhang PetscInt nis = iscoloring->n,*rowhit,*columnsforrow,bs = 1; 745b1683b59SHong Zhang IS *isa; 746b1683b59SHong Zhang PetscBool done; 747b1683b59SHong Zhang PetscBool flg1,flg2; 748b1683b59SHong Zhang 749b1683b59SHong Zhang PetscFunctionBegin; 750b1683b59SHong Zhang if (!mat->assembled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix must be assembled by calls to MatAssemblyBegin/End();"); 751b1683b59SHong Zhang 752b1683b59SHong Zhang ierr = ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);CHKERRQ(ierr); 753b1683b59SHong Zhang /* this is ugly way to get blocksize but cannot call MatGetBlockSize() because AIJ can have bs > 1 */ 754b1683b59SHong Zhang ierr = PetscTypeCompare((PetscObject)mat,MATSEQBAIJ,&flg1);CHKERRQ(ierr); 755b1683b59SHong Zhang ierr = PetscTypeCompare((PetscObject)mat,MATMPIBAIJ,&flg2);CHKERRQ(ierr); 756b1683b59SHong Zhang if (flg1 || flg2) { 757b1683b59SHong Zhang ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr); 758b1683b59SHong Zhang } 759b1683b59SHong Zhang 760b1683b59SHong Zhang N = mat->cmap->N/bs; 761b1683b59SHong Zhang c->M = mat->rmap->N/bs; /* set total rows, columns and local rows */ 762b1683b59SHong Zhang c->N = mat->cmap->N/bs; 763b1683b59SHong Zhang c->m = mat->rmap->N/bs; 764b1683b59SHong Zhang c->rstart = 0; 765b1683b59SHong Zhang 766b1683b59SHong Zhang c->ncolors = nis; 767b1683b59SHong Zhang ierr = PetscMalloc(nis*sizeof(PetscInt),&c->ncolumns);CHKERRQ(ierr); 768b1683b59SHong Zhang ierr = PetscMalloc(nis*sizeof(PetscInt*),&c->columns);CHKERRQ(ierr); 769b1683b59SHong Zhang ierr = PetscMalloc(c->m*sizeof(PetscInt),&c->nrows);CHKERRQ(ierr); 770b1683b59SHong Zhang ierr = PetscMalloc(c->m*sizeof(PetscInt*),&c->rows);CHKERRQ(ierr); 771b1683b59SHong Zhang ierr = PetscMalloc(c->m*sizeof(PetscInt*),&c->columnsforrow);CHKERRQ(ierr); 772b1683b59SHong Zhang 773b1683b59SHong Zhang ierr = MatGetColumnIJ(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&done);CHKERRQ(ierr); 774b1683b59SHong Zhang if (!done) SETERRQ1(((PetscObject)mat)->comm,PETSC_ERR_SUP,"MatGetColumnIJ() not supported for matrix type %s",((PetscObject)mat)->type_name); 775b1683b59SHong Zhang 776b1683b59SHong Zhang ierr = PetscMalloc((c->m+1)*sizeof(PetscInt),&rowhit);CHKERRQ(ierr); 777b1683b59SHong Zhang ierr = PetscMalloc((c->m+1)*sizeof(PetscInt),&columnsforrow);CHKERRQ(ierr); 778b1683b59SHong Zhang 779b1683b59SHong Zhang for (i=0; i<nis; i++) { 780b1683b59SHong Zhang ierr = ISGetLocalSize(isa[i],&n);CHKERRQ(ierr); 781b1683b59SHong Zhang ierr = ISGetIndices(isa[i],&is);CHKERRQ(ierr); 782b1683b59SHong Zhang c->ncolumns[i] = n; 783b1683b59SHong Zhang if (n) { 784b1683b59SHong Zhang ierr = PetscMalloc(n*sizeof(PetscInt),&c->columns[i]);CHKERRQ(ierr); 785b1683b59SHong Zhang ierr = PetscMemcpy(c->columns[i],is,n*sizeof(PetscInt));CHKERRQ(ierr); 786b1683b59SHong Zhang } else { 787b1683b59SHong Zhang c->columns[i] = 0; 788b1683b59SHong Zhang } 789b1683b59SHong Zhang 790b1683b59SHong Zhang /* fast, crude version requires O(N*N) work */ 791b1683b59SHong Zhang ierr = PetscMemzero(rowhit,c->m*sizeof(PetscInt));CHKERRQ(ierr); 792b1683b59SHong Zhang /* loop over columns*/ 793b1683b59SHong Zhang for (j=0; j<n; j++) { 794b1683b59SHong Zhang col = is[j]; 795b1683b59SHong Zhang rows = cj + ci[col]; 796b1683b59SHong Zhang m = ci[col+1] - ci[col]; 797b1683b59SHong Zhang /* loop over columns marking them in rowhit */ 798b1683b59SHong Zhang for (k=0; k<m; k++) { 799b1683b59SHong Zhang rowhit[*rows++] = col + 1; 800b1683b59SHong Zhang } 801b1683b59SHong Zhang } 802b1683b59SHong Zhang /* count the number of hits */ 803b1683b59SHong Zhang nrows = 0; 804b1683b59SHong Zhang for (j=0; j<c->m; j++) { 805b1683b59SHong Zhang if (rowhit[j]) nrows++; 806b1683b59SHong Zhang } 807b1683b59SHong Zhang c->nrows[i] = nrows; 808b1683b59SHong Zhang ierr = PetscMalloc((nrows+1)*sizeof(PetscInt),&c->rows[i]);CHKERRQ(ierr); 809b1683b59SHong Zhang ierr = PetscMalloc((nrows+1)*sizeof(PetscInt),&c->columnsforrow[i]);CHKERRQ(ierr); 810b1683b59SHong Zhang nrows = 0; 811b1683b59SHong Zhang for (j=0; j<c->m; j++) { 812b1683b59SHong Zhang if (rowhit[j]) { 813b1683b59SHong Zhang c->rows[i][nrows] = j; 814b1683b59SHong Zhang c->columnsforrow[i][nrows] = rowhit[j] - 1; 815b1683b59SHong Zhang nrows++; 816b1683b59SHong Zhang } 817b1683b59SHong Zhang } 818b1683b59SHong Zhang ierr = ISRestoreIndices(isa[i],&is);CHKERRQ(ierr); 819b1683b59SHong Zhang } 820b1683b59SHong Zhang ierr = MatRestoreColumnIJ(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&done);CHKERRQ(ierr); 821b1683b59SHong Zhang 822b1683b59SHong Zhang ierr = PetscFree(rowhit);CHKERRQ(ierr); 823b1683b59SHong Zhang ierr = PetscFree(columnsforrow);CHKERRQ(ierr); 824b1683b59SHong Zhang ierr = ISColoringRestoreIS(iscoloring,&isa);CHKERRQ(ierr); 825b1683b59SHong Zhang PetscFunctionReturn(0); 826b1683b59SHong Zhang } 827