1be1d678aSKris Buschelman #define PETSCMAT_DLL 2be1d678aSKris Buschelman 3d50806bdSBarry Smith /* 42499ec78SHong Zhang Defines matrix-matrix product routines for pairs of SeqAIJ matrices 5d50806bdSBarry Smith C = A * B 6d50806bdSBarry Smith */ 7d50806bdSBarry Smith 8c1f4806aSKris Buschelman #include "src/mat/impls/aij/seq/aij.h" /*I "petscmat.h" I*/ 970f19b1fSKris Buschelman #include "src/mat/utils/freespace.h" 10be0fcf8dSHong Zhang #include "petscbt.h" 114ae313f4SHong Zhang #include "src/mat/impls/dense/seq/dense.h" /*I "petscmat.h" I*/ 126284ec50SHong Zhang 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 } 265c66b693SKris Buschelman 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; 36899cda47SBarry 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; 4058c24d83SHong Zhang 4158c24d83SHong Zhang PetscFunctionBegin; 4258c24d83SHong Zhang /* Set up */ 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) { 75a1a86e44SBarry Smith ierr = PetscFreeSpaceGet(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 */ 10058c24d83SHong Zhang ierr = MatCreateSeqAIJWithArrays(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); 105*e6b907acSBarry Smith c->free_a = PETSC_TRUE; 106*e6b907acSBarry Smith c->free_ij = PETSC_TRUE; 10758c24d83SHong Zhang c->nonew = 0; 10858c24d83SHong Zhang 109be0fcf8dSHong Zhang if (nspacedouble){ 110ae15b995SBarry Smith ierr = PetscInfo5((*C),"nspacedouble:%D, nnz(A):%D, nnz(B):%D, fill:%G, nnz(C):%D\n",nspacedouble,ai[am],bi[bm],fill,ci[am]);CHKERRQ(ierr); 111be0fcf8dSHong Zhang } 11258c24d83SHong Zhang PetscFunctionReturn(0); 11358c24d83SHong Zhang } 114d50806bdSBarry Smith 1155c66b693SKris Buschelman 11626be0446SHong Zhang #undef __FUNCT__ 11726be0446SHong Zhang #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqAIJ" 118dfbe8321SBarry Smith PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C) 119d50806bdSBarry Smith { 120dfbe8321SBarry Smith PetscErrorCode ierr; 12138baddfdSBarry Smith PetscInt flops=0; 122d50806bdSBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 123d50806bdSBarry Smith Mat_SeqAIJ *b = (Mat_SeqAIJ *)B->data; 124d50806bdSBarry Smith Mat_SeqAIJ *c = (Mat_SeqAIJ *)C->data; 12538baddfdSBarry Smith PetscInt *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j; 126899cda47SBarry Smith PetscInt am=A->rmap.N,cm=C->rmap.N; 127c124e916SHong Zhang PetscInt i,j,k,anzi,bnzi,cnzi,brow,nextb; 128c124e916SHong Zhang MatScalar *aa=a->a,*ba=b->a,*baj,*ca=c->a; 129d50806bdSBarry Smith 130d50806bdSBarry Smith PetscFunctionBegin; 131c124e916SHong Zhang /* clean old values in C */ 132c124e916SHong Zhang ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr); 133d50806bdSBarry Smith /* Traverse A row-wise. */ 134d50806bdSBarry Smith /* Build the ith row in C by summing over nonzero columns in A, */ 135d50806bdSBarry Smith /* the rows of B corresponding to nonzeros of A. */ 136d50806bdSBarry Smith for (i=0;i<am;i++) { 137d50806bdSBarry Smith anzi = ai[i+1] - ai[i]; 138d50806bdSBarry Smith for (j=0;j<anzi;j++) { 139d50806bdSBarry Smith brow = *aj++; 140d50806bdSBarry Smith bnzi = bi[brow+1] - bi[brow]; 141d50806bdSBarry Smith bjj = bj + bi[brow]; 142d50806bdSBarry Smith baj = ba + bi[brow]; 143c124e916SHong Zhang nextb = 0; 144c124e916SHong Zhang for (k=0; nextb<bnzi; k++) { 145c124e916SHong Zhang if (cj[k] == bjj[nextb]){ /* ccol == bcol */ 146c124e916SHong Zhang ca[k] += (*aa)*baj[nextb++]; 147c124e916SHong Zhang } 148d50806bdSBarry Smith } 149d50806bdSBarry Smith flops += 2*bnzi; 150d50806bdSBarry Smith aa++; 151d50806bdSBarry Smith } 152d50806bdSBarry Smith cnzi = ci[i+1] - ci[i]; 153d50806bdSBarry Smith ca += cnzi; 154d50806bdSBarry Smith cj += cnzi; 155d50806bdSBarry Smith } 156716bacf3SKris Buschelman ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 157716bacf3SKris Buschelman ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 158716bacf3SKris Buschelman 159d50806bdSBarry Smith ierr = PetscLogFlops(flops);CHKERRQ(ierr); 160d50806bdSBarry Smith PetscFunctionReturn(0); 161d50806bdSBarry Smith } 162bc011b1eSHong Zhang 163bc011b1eSHong Zhang 164bc011b1eSHong Zhang #undef __FUNCT__ 165bc011b1eSHong Zhang #define __FUNCT__ "MatMatMultTranspose_SeqAIJ_SeqAIJ" 166bc011b1eSHong Zhang PetscErrorCode MatMatMultTranspose_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) { 167bc011b1eSHong Zhang PetscErrorCode ierr; 168bc011b1eSHong Zhang 169bc011b1eSHong Zhang PetscFunctionBegin; 170bc011b1eSHong Zhang if (scall == MAT_INITIAL_MATRIX){ 171bc011b1eSHong Zhang ierr = MatMatMultTransposeSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr); 172bc011b1eSHong Zhang } 173bc011b1eSHong Zhang ierr = MatMatMultTransposeNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr); 174bc011b1eSHong Zhang PetscFunctionReturn(0); 175bc011b1eSHong Zhang } 176bc011b1eSHong Zhang 177bc011b1eSHong Zhang #undef __FUNCT__ 178bc011b1eSHong Zhang #define __FUNCT__ "MatMatMultTransposeSymbolic_SeqAIJ_SeqAIJ" 179bc011b1eSHong Zhang PetscErrorCode MatMatMultTransposeSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C) 180bc011b1eSHong Zhang { 181bc011b1eSHong Zhang PetscErrorCode ierr; 182bc011b1eSHong Zhang Mat At; 18338baddfdSBarry Smith PetscInt *ati,*atj; 184bc011b1eSHong Zhang 185bc011b1eSHong Zhang PetscFunctionBegin; 186bc011b1eSHong Zhang /* create symbolic At */ 187bc011b1eSHong Zhang ierr = MatGetSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr); 188899cda47SBarry Smith ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,A->cmap.n,A->rmap.n,ati,atj,PETSC_NULL,&At);CHKERRQ(ierr); 189bc011b1eSHong Zhang 190bc011b1eSHong Zhang /* get symbolic C=At*B */ 191bc011b1eSHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(At,B,fill,C);CHKERRQ(ierr); 192bc011b1eSHong Zhang 193bc011b1eSHong Zhang /* clean up */ 194bc011b1eSHong Zhang ierr = MatDestroy(At);CHKERRQ(ierr); 195bc011b1eSHong Zhang ierr = MatRestoreSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr); 196bc011b1eSHong Zhang 197bc011b1eSHong Zhang PetscFunctionReturn(0); 198bc011b1eSHong Zhang } 199bc011b1eSHong Zhang 200bc011b1eSHong Zhang #undef __FUNCT__ 201bc011b1eSHong Zhang #define __FUNCT__ "MatMatMultTransposeNumeric_SeqAIJ_SeqAIJ" 202bc011b1eSHong Zhang PetscErrorCode MatMatMultTransposeNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C) 203bc011b1eSHong Zhang { 204bc011b1eSHong Zhang PetscErrorCode ierr; 2050fbc74f4SHong Zhang Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data; 206899cda47SBarry Smith PetscInt am=A->rmap.n,anzi,*ai=a->i,*aj=a->j,*bi=b->i,*bj,bnzi,nextb; 207899cda47SBarry Smith PetscInt cm=C->rmap.n,*ci=c->i,*cj=c->j,crow,*cjj,i,j,k,flops=0; 2080fbc74f4SHong Zhang MatScalar *aa=a->a,*ba,*ca=c->a,*caj; 209bc011b1eSHong Zhang 210bc011b1eSHong Zhang PetscFunctionBegin; 211bc011b1eSHong Zhang /* clear old values in C */ 212bc011b1eSHong Zhang ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr); 213bc011b1eSHong Zhang 214bc011b1eSHong Zhang /* compute A^T*B using outer product (A^T)[:,i]*B[i,:] */ 215bc011b1eSHong Zhang for (i=0;i<am;i++) { 216bc011b1eSHong Zhang bj = b->j + bi[i]; 217bc011b1eSHong Zhang ba = b->a + bi[i]; 218bc011b1eSHong Zhang bnzi = bi[i+1] - bi[i]; 219bc011b1eSHong Zhang anzi = ai[i+1] - ai[i]; 220bc011b1eSHong Zhang for (j=0; j<anzi; j++) { 221bc011b1eSHong Zhang nextb = 0; 2220fbc74f4SHong Zhang crow = *aj++; 223bc011b1eSHong Zhang cjj = cj + ci[crow]; 224bc011b1eSHong Zhang caj = ca + ci[crow]; 225bc011b1eSHong Zhang /* perform sparse axpy operation. Note cjj includes bj. */ 226bc011b1eSHong Zhang for (k=0; nextb<bnzi; k++) { 2270fbc74f4SHong Zhang if (cjj[k] == *(bj+nextb)) { /* ccol == bcol */ 2280fbc74f4SHong Zhang caj[k] += (*aa)*(*(ba+nextb)); 229bc011b1eSHong Zhang nextb++; 230bc011b1eSHong Zhang } 231bc011b1eSHong Zhang } 232bc011b1eSHong Zhang flops += 2*bnzi; 2330fbc74f4SHong Zhang aa++; 234bc011b1eSHong Zhang } 235bc011b1eSHong Zhang } 236bc011b1eSHong Zhang 237bc011b1eSHong Zhang /* Assemble the final matrix and clean up */ 238bc011b1eSHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 239bc011b1eSHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 240bc011b1eSHong Zhang ierr = PetscLogFlops(flops);CHKERRQ(ierr); 241bc011b1eSHong Zhang PetscFunctionReturn(0); 242bc011b1eSHong Zhang } 2439123193aSHong Zhang 2449123193aSHong Zhang #undef __FUNCT__ 2459123193aSHong Zhang #define __FUNCT__ "MatMatMult_SeqAIJ_SeqDense" 2469123193aSHong Zhang PetscErrorCode MatMatMult_SeqAIJ_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 2479123193aSHong Zhang { 2489123193aSHong Zhang PetscErrorCode ierr; 2499123193aSHong Zhang 2509123193aSHong Zhang PetscFunctionBegin; 2519123193aSHong Zhang if (scall == MAT_INITIAL_MATRIX){ 2529123193aSHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqDense(A,B,fill,C);CHKERRQ(ierr); 2539123193aSHong Zhang } 2549123193aSHong Zhang ierr = MatMatMultNumeric_SeqAIJ_SeqDense(A,B,*C);CHKERRQ(ierr); 2559123193aSHong Zhang PetscFunctionReturn(0); 2569123193aSHong Zhang } 257edd81eecSMatthew Knepley 2589123193aSHong Zhang #undef __FUNCT__ 2599123193aSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqDense" 2609123193aSHong Zhang PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C) 2619123193aSHong Zhang { 2629123193aSHong Zhang PetscErrorCode ierr; 2639123193aSHong Zhang 2649123193aSHong Zhang PetscFunctionBegin; 2659123193aSHong Zhang ierr = MatMatMultSymbolic_SeqDense_SeqDense(A,B,0.0,C); 2669123193aSHong Zhang PetscFunctionReturn(0); 2679123193aSHong Zhang } 2689123193aSHong Zhang 2699123193aSHong Zhang #undef __FUNCT__ 2709123193aSHong Zhang #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqDense" 2719123193aSHong Zhang PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqDense(Mat A,Mat B,Mat C) 2729123193aSHong Zhang { 273f32d5d43SBarry Smith Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data; 2749123193aSHong Zhang PetscErrorCode ierr; 275f32d5d43SBarry Smith PetscScalar *b,*c,r1,r2,r3,r4,*aa,*b1,*b2,*b3,*b4; 276f32d5d43SBarry Smith PetscInt cm=C->rmap.n, cn=B->cmap.n, bm=B->rmap.n, col, i,j,n,*aj, am = A->rmap.n; 277f32d5d43SBarry Smith PetscInt am2 = 2*am, am3 = 3*am, bm4 = 4*bm,colam; 2789123193aSHong Zhang 2799123193aSHong Zhang PetscFunctionBegin; 280f32d5d43SBarry Smith if (!cm || !cn) PetscFunctionReturn(0); 281f32d5d43SBarry Smith if (bm != A->cmap.n) SETERRQ2(PETSC_ERR_ARG_SIZ,"Number columns in A %D not equal rows in B %D\n",A->cmap.n,bm); 282f32d5d43SBarry Smith if (A->rmap.n != C->rmap.n) SETERRQ2(PETSC_ERR_ARG_SIZ,"Number rows in C %D not equal rows in A %D\n",C->rmap.n,A->rmap.n); 283f32d5d43SBarry Smith if (B->cmap.n != C->cmap.n) SETERRQ2(PETSC_ERR_ARG_SIZ,"Number columns in B %D not equal columns in C %D\n",B->cmap.n,C->cmap.n); 284f32d5d43SBarry Smith ierr = MatGetArray(B,&b);CHKERRQ(ierr); 285f32d5d43SBarry Smith ierr = MatGetArray(C,&c);CHKERRQ(ierr); 286f32d5d43SBarry Smith b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm; 287f32d5d43SBarry Smith for (col=0; col<cn-4; col += 4){ /* over columns of C */ 288f32d5d43SBarry Smith colam = col*am; 289f32d5d43SBarry Smith for (i=0; i<am; i++) { /* over rows of C in those columns */ 290f32d5d43SBarry Smith r1 = r2 = r3 = r4 = 0.0; 291f32d5d43SBarry Smith n = a->i[i+1] - a->i[i]; 292f32d5d43SBarry Smith aj = a->j + a->i[i]; 293f32d5d43SBarry Smith aa = a->a + a->i[i]; 294f32d5d43SBarry Smith for (j=0; j<n; j++) { 295f32d5d43SBarry Smith r1 += (*aa)*b1[*aj]; 296f32d5d43SBarry Smith r2 += (*aa)*b2[*aj]; 297f32d5d43SBarry Smith r3 += (*aa)*b3[*aj]; 298f32d5d43SBarry Smith r4 += (*aa++)*b4[*aj++]; 2999123193aSHong Zhang } 300f32d5d43SBarry Smith c[colam + i] = r1; 301f32d5d43SBarry Smith c[colam + am + i] = r2; 302f32d5d43SBarry Smith c[colam + am2 + i] = r3; 303f32d5d43SBarry Smith c[colam + am3 + i] = r4; 304f32d5d43SBarry Smith } 305f32d5d43SBarry Smith b1 += bm4; 306f32d5d43SBarry Smith b2 += bm4; 307f32d5d43SBarry Smith b3 += bm4; 308f32d5d43SBarry Smith b4 += bm4; 309f32d5d43SBarry Smith } 310f32d5d43SBarry Smith for (;col<cn; col++){ /* over extra columns of C */ 311f32d5d43SBarry Smith for (i=0; i<am; i++) { /* over rows of C in those columns */ 312f32d5d43SBarry Smith r1 = 0.0; 313f32d5d43SBarry Smith n = a->i[i+1] - a->i[i]; 314f32d5d43SBarry Smith aj = a->j + a->i[i]; 315f32d5d43SBarry Smith aa = a->a + a->i[i]; 316f32d5d43SBarry Smith 317f32d5d43SBarry Smith for (j=0; j<n; j++) { 318f32d5d43SBarry Smith r1 += (*aa++)*b1[*aj++]; 319f32d5d43SBarry Smith } 320f32d5d43SBarry Smith c[col*am + i] = r1; 321f32d5d43SBarry Smith } 322f32d5d43SBarry Smith b1 += bm; 323f32d5d43SBarry Smith } 3244324174eSBarry Smith ierr = PetscLogFlops(cn*(2*a->nz - A->rmap.n));CHKERRQ(ierr); 325f32d5d43SBarry Smith ierr = MatRestoreArray(B,&b);CHKERRQ(ierr); 326f32d5d43SBarry Smith ierr = MatRestoreArray(C,&c);CHKERRQ(ierr); 327f32d5d43SBarry Smith ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 328f32d5d43SBarry Smith ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 329f32d5d43SBarry Smith PetscFunctionReturn(0); 330f32d5d43SBarry Smith } 331f32d5d43SBarry Smith 332f32d5d43SBarry Smith /* 3334324174eSBarry Smith Note very similar to MatMult_SeqAIJ(), should generate both codes from same base 334f32d5d43SBarry Smith */ 335f32d5d43SBarry Smith #undef __FUNCT__ 336f32d5d43SBarry Smith #define __FUNCT__ "MatMatMultNumericAdd_SeqAIJ_SeqDense" 337f32d5d43SBarry Smith PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat A,Mat B,Mat C) 338f32d5d43SBarry Smith { 339f32d5d43SBarry Smith Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data; 340f32d5d43SBarry Smith PetscErrorCode ierr; 341f32d5d43SBarry Smith PetscScalar *b,*c,r1,r2,r3,r4,*aa,*b1,*b2,*b3,*b4; 3424324174eSBarry 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; 3434324174eSBarry Smith PetscInt am2 = 2*am, am3 = 3*am, bm4 = 4*bm,colam,*ridx; 344f32d5d43SBarry Smith 345f32d5d43SBarry Smith PetscFunctionBegin; 346f32d5d43SBarry Smith if (!cm || !cn) PetscFunctionReturn(0); 347f32d5d43SBarry Smith ierr = MatGetArray(B,&b);CHKERRQ(ierr); 348f32d5d43SBarry Smith ierr = MatGetArray(C,&c);CHKERRQ(ierr); 349f32d5d43SBarry Smith b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm; 3504324174eSBarry Smith 3514324174eSBarry Smith if (a->compressedrow.use){ /* use compressed row format */ 3524324174eSBarry Smith for (col=0; col<cn-4; col += 4){ /* over columns of C */ 3534324174eSBarry Smith colam = col*am; 3544324174eSBarry Smith arm = a->compressedrow.nrows; 3554324174eSBarry Smith ii = a->compressedrow.i; 3564324174eSBarry Smith ridx = a->compressedrow.rindex; 3574324174eSBarry Smith for (i=0; i<arm; i++) { /* over rows of C in those columns */ 3584324174eSBarry Smith r1 = r2 = r3 = r4 = 0.0; 3594324174eSBarry Smith n = ii[i+1] - ii[i]; 3604324174eSBarry Smith aj = a->j + ii[i]; 3614324174eSBarry Smith aa = a->a + ii[i]; 3624324174eSBarry Smith for (j=0; j<n; j++) { 3634324174eSBarry Smith r1 += (*aa)*b1[*aj]; 3644324174eSBarry Smith r2 += (*aa)*b2[*aj]; 3654324174eSBarry Smith r3 += (*aa)*b3[*aj]; 3664324174eSBarry Smith r4 += (*aa++)*b4[*aj++]; 3674324174eSBarry Smith } 3684324174eSBarry Smith c[colam + ridx[i]] += r1; 3694324174eSBarry Smith c[colam + am + ridx[i]] += r2; 3704324174eSBarry Smith c[colam + am2 + ridx[i]] += r3; 3714324174eSBarry Smith c[colam + am3 + ridx[i]] += r4; 3724324174eSBarry Smith } 3734324174eSBarry Smith b1 += bm4; 3744324174eSBarry Smith b2 += bm4; 3754324174eSBarry Smith b3 += bm4; 3764324174eSBarry Smith b4 += bm4; 3774324174eSBarry Smith } 3784324174eSBarry Smith for (;col<cn; col++){ /* over extra columns of C */ 3794324174eSBarry Smith colam = col*am; 3804324174eSBarry Smith arm = a->compressedrow.nrows; 3814324174eSBarry Smith ii = a->compressedrow.i; 3824324174eSBarry Smith ridx = a->compressedrow.rindex; 3834324174eSBarry Smith for (i=0; i<arm; i++) { /* over rows of C in those columns */ 3844324174eSBarry Smith r1 = 0.0; 3854324174eSBarry Smith n = ii[i+1] - ii[i]; 3864324174eSBarry Smith aj = a->j + ii[i]; 3874324174eSBarry Smith aa = a->a + ii[i]; 3884324174eSBarry Smith 3894324174eSBarry Smith for (j=0; j<n; j++) { 3904324174eSBarry Smith r1 += (*aa++)*b1[*aj++]; 3914324174eSBarry Smith } 3924324174eSBarry Smith c[col*am + ridx[i]] += r1; 3934324174eSBarry Smith } 3944324174eSBarry Smith b1 += bm; 3954324174eSBarry Smith } 3964324174eSBarry Smith } else { 397f32d5d43SBarry Smith for (col=0; col<cn-4; col += 4){ /* over columns of C */ 398f32d5d43SBarry Smith colam = col*am; 399f32d5d43SBarry Smith for (i=0; i<am; i++) { /* over rows of C in those columns */ 400f32d5d43SBarry Smith r1 = r2 = r3 = r4 = 0.0; 401f32d5d43SBarry Smith n = a->i[i+1] - a->i[i]; 402f32d5d43SBarry Smith aj = a->j + a->i[i]; 403f32d5d43SBarry Smith aa = a->a + a->i[i]; 404f32d5d43SBarry Smith for (j=0; j<n; j++) { 405f32d5d43SBarry Smith r1 += (*aa)*b1[*aj]; 406f32d5d43SBarry Smith r2 += (*aa)*b2[*aj]; 407f32d5d43SBarry Smith r3 += (*aa)*b3[*aj]; 408f32d5d43SBarry Smith r4 += (*aa++)*b4[*aj++]; 409f32d5d43SBarry Smith } 410f32d5d43SBarry Smith c[colam + i] += r1; 411f32d5d43SBarry Smith c[colam + am + i] += r2; 412f32d5d43SBarry Smith c[colam + am2 + i] += r3; 413f32d5d43SBarry Smith c[colam + am3 + i] += r4; 414f32d5d43SBarry Smith } 415f32d5d43SBarry Smith b1 += bm4; 416f32d5d43SBarry Smith b2 += bm4; 417f32d5d43SBarry Smith b3 += bm4; 418f32d5d43SBarry Smith b4 += bm4; 419f32d5d43SBarry Smith } 420f32d5d43SBarry Smith for (;col<cn; col++){ /* over extra columns of C */ 421f32d5d43SBarry Smith for (i=0; i<am; i++) { /* over rows of C in those columns */ 422f32d5d43SBarry Smith r1 = 0.0; 423f32d5d43SBarry Smith n = a->i[i+1] - a->i[i]; 424f32d5d43SBarry Smith aj = a->j + a->i[i]; 425f32d5d43SBarry Smith aa = a->a + a->i[i]; 426f32d5d43SBarry Smith 427f32d5d43SBarry Smith for (j=0; j<n; j++) { 428f32d5d43SBarry Smith r1 += (*aa++)*b1[*aj++]; 429f32d5d43SBarry Smith } 430f32d5d43SBarry Smith c[col*am + i] += r1; 431f32d5d43SBarry Smith } 432f32d5d43SBarry Smith b1 += bm; 433f32d5d43SBarry Smith } 4344324174eSBarry Smith } 4354324174eSBarry Smith ierr = PetscLogFlops(cn*2*a->nz);CHKERRQ(ierr); 436f32d5d43SBarry Smith ierr = MatRestoreArray(B,&b);CHKERRQ(ierr); 437f32d5d43SBarry Smith ierr = MatRestoreArray(C,&c);CHKERRQ(ierr); 4389123193aSHong Zhang PetscFunctionReturn(0); 4399123193aSHong Zhang } 440