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 177d2853540SHong Zhang /* This routine is not used. Should be removed! */ 178bc011b1eSHong Zhang #undef __FUNCT__ 1796fc122caSHong Zhang #define __FUNCT__ "MatMatTransposeMult_SeqAIJ_SeqAIJ" 1806fc122caSHong Zhang PetscErrorCode MatMatTransposeMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 1815df89d91SHong Zhang { 182bc011b1eSHong Zhang PetscErrorCode ierr; 183bc011b1eSHong Zhang 184bc011b1eSHong Zhang PetscFunctionBegin; 185bc011b1eSHong Zhang if (scall == MAT_INITIAL_MATRIX){ 1866fc122caSHong Zhang ierr = MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr); 187bc011b1eSHong Zhang } 1886fc122caSHong Zhang ierr = MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr); 189bc011b1eSHong Zhang PetscFunctionReturn(0); 190bc011b1eSHong Zhang } 191bc011b1eSHong Zhang 192bc011b1eSHong Zhang #undef __FUNCT__ 1932128a86cSHong Zhang #define __FUNCT__ "PetscContainerDestroy_Mat_MatMatMultTrans" 1942128a86cSHong Zhang PetscErrorCode PetscContainerDestroy_Mat_MatMatMultTrans(void *ptr) 1952128a86cSHong Zhang { 1962128a86cSHong Zhang PetscErrorCode ierr; 1972128a86cSHong Zhang Mat_MatMatMultTrans *multtrans=(Mat_MatMatMultTrans*)ptr; 1982128a86cSHong Zhang 1992128a86cSHong Zhang PetscFunctionBegin; 200b9af6bddSHong Zhang ierr = MatTransposeColoringDestroy(&multtrans->matcoloring);CHKERRQ(ierr); 2012128a86cSHong Zhang ierr = MatDestroy(&multtrans->Bt_den);CHKERRQ(ierr); 2022128a86cSHong Zhang ierr = MatDestroy(&multtrans->ABt_den);CHKERRQ(ierr); 2032128a86cSHong Zhang ierr = PetscFree(multtrans);CHKERRQ(ierr); 2042128a86cSHong Zhang PetscFunctionReturn(0); 2052128a86cSHong Zhang } 2062128a86cSHong Zhang 2072128a86cSHong Zhang #undef __FUNCT__ 2082128a86cSHong Zhang #define __FUNCT__ "MatDestroy_SeqAIJ_MatMatMultTrans" 2092128a86cSHong Zhang PetscErrorCode MatDestroy_SeqAIJ_MatMatMultTrans(Mat A) 2102128a86cSHong Zhang { 2112128a86cSHong Zhang PetscErrorCode ierr; 2122128a86cSHong Zhang PetscContainer container; 2132128a86cSHong Zhang Mat_MatMatMultTrans *multtrans=PETSC_NULL; 2142128a86cSHong Zhang 2152128a86cSHong Zhang PetscFunctionBegin; 2162128a86cSHong Zhang ierr = PetscObjectQuery((PetscObject)A,"Mat_MatMatMultTrans",(PetscObject *)&container);CHKERRQ(ierr); 2172128a86cSHong Zhang if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit"); 2182128a86cSHong Zhang ierr = PetscContainerGetPointer(container,(void **)&multtrans);CHKERRQ(ierr); 2192128a86cSHong Zhang A->ops->destroy = multtrans->destroy; 2202128a86cSHong Zhang if (A->ops->destroy) { 2212128a86cSHong Zhang ierr = (*A->ops->destroy)(A);CHKERRQ(ierr); 2222128a86cSHong Zhang } 2232128a86cSHong Zhang ierr = PetscObjectCompose((PetscObject)A,"Mat_MatMatMultTrans",0);CHKERRQ(ierr); 2242128a86cSHong Zhang PetscFunctionReturn(0); 2252128a86cSHong Zhang } 2262128a86cSHong Zhang 2272128a86cSHong Zhang #undef __FUNCT__ 2286fc122caSHong Zhang #define __FUNCT__ "MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ" 2296fc122caSHong Zhang PetscErrorCode MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C) 230bc011b1eSHong Zhang { 2315df89d91SHong Zhang PetscErrorCode ierr; 23281d82fe4SHong Zhang Mat Bt; 23381d82fe4SHong Zhang PetscInt *bti,*btj; 2342128a86cSHong Zhang Mat_MatMatMultTrans *multtrans; 2352128a86cSHong Zhang PetscContainer container; 236d2853540SHong Zhang PetscLogDouble t0,tf,etime2=0.0; 237d2853540SHong Zhang 23881d82fe4SHong Zhang PetscFunctionBegin; 239d2853540SHong Zhang ierr = PetscGetTime(&t0);CHKERRQ(ierr); 24081d82fe4SHong Zhang /* create symbolic Bt */ 24181d82fe4SHong Zhang ierr = MatGetSymbolicTranspose_SeqAIJ(B,&bti,&btj);CHKERRQ(ierr); 24281d82fe4SHong Zhang ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,B->cmap->n,B->rmap->n,bti,btj,PETSC_NULL,&Bt);CHKERRQ(ierr); 24381d82fe4SHong Zhang 24481d82fe4SHong Zhang /* get symbolic C=A*Bt */ 24581d82fe4SHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,Bt,fill,C);CHKERRQ(ierr); 24681d82fe4SHong Zhang 2472128a86cSHong Zhang /* create a supporting struct for reuse intermidiate dense matrices with matcoloring */ 2482128a86cSHong Zhang ierr = PetscNew(Mat_MatMatMultTrans,&multtrans);CHKERRQ(ierr); 2492128a86cSHong Zhang 2502128a86cSHong Zhang /* attach the supporting struct to C */ 2512128a86cSHong Zhang ierr = PetscContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr); 2522128a86cSHong Zhang ierr = PetscContainerSetPointer(container,multtrans);CHKERRQ(ierr); 2532128a86cSHong Zhang ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_Mat_MatMatMultTrans);CHKERRQ(ierr); 2542128a86cSHong Zhang ierr = PetscObjectCompose((PetscObject)(*C),"Mat_MatMatMultTrans",(PetscObject)container);CHKERRQ(ierr); 2552128a86cSHong Zhang ierr = PetscContainerDestroy(&container);CHKERRQ(ierr); 2562128a86cSHong Zhang 2572128a86cSHong Zhang multtrans->usecoloring = PETSC_FALSE; 2582128a86cSHong Zhang multtrans->destroy = (*C)->ops->destroy; 2592128a86cSHong Zhang (*C)->ops->destroy = MatDestroy_SeqAIJ_MatMatMultTrans; 2602128a86cSHong Zhang 261d2853540SHong Zhang ierr = PetscGetTime(&tf);CHKERRQ(ierr); 262d2853540SHong Zhang etime2 += tf - t0; 263d2853540SHong Zhang 264f400d4cbSHong Zhang ierr = PetscOptionsGetBool(PETSC_NULL,"-matmattransmult_color",&multtrans->usecoloring,PETSC_NULL);CHKERRQ(ierr); 2652128a86cSHong Zhang if (multtrans->usecoloring){ 266b9af6bddSHong Zhang /* Create MatTransposeColoring from symbolic C=A*B^T */ 267b9af6bddSHong Zhang MatTransposeColoring matcoloring; 2682128a86cSHong Zhang ISColoring iscoloring; 2692128a86cSHong Zhang Mat Bt_dense,C_dense; 270d2853540SHong Zhang PetscLogDouble etime0=0.0,etime01=0.0,etime1=0.0; 271e8354b3eSHong Zhang 272e8354b3eSHong Zhang ierr = PetscGetTime(&t0);CHKERRQ(ierr); 273502de53fSHong Zhang ierr = MatGetColoring(*C,MATCOLORINGLF,&iscoloring);CHKERRQ(ierr); 274d2853540SHong Zhang ierr = PetscGetTime(&tf);CHKERRQ(ierr); 275d2853540SHong Zhang etime0 += tf - t0; 276d2853540SHong Zhang 277d2853540SHong Zhang ierr = PetscGetTime(&t0);CHKERRQ(ierr); 278b9af6bddSHong Zhang ierr = MatTransposeColoringCreate(*C,iscoloring,&matcoloring);CHKERRQ(ierr); 2792128a86cSHong Zhang multtrans->matcoloring = matcoloring; 2802128a86cSHong Zhang ierr = ISColoringDestroy(&iscoloring);CHKERRQ(ierr); 281e8354b3eSHong Zhang ierr = PetscGetTime(&tf);CHKERRQ(ierr); 282d2853540SHong Zhang etime01 += tf - t0; 2832128a86cSHong Zhang 284e8354b3eSHong Zhang ierr = PetscGetTime(&t0);CHKERRQ(ierr); 2852128a86cSHong Zhang /* Create Bt_dense and C_dense = A*Bt_dense */ 2862128a86cSHong Zhang ierr = MatCreate(PETSC_COMM_SELF,&Bt_dense);CHKERRQ(ierr); 2872128a86cSHong Zhang ierr = MatSetSizes(Bt_dense,A->cmap->n,matcoloring->ncolors,A->cmap->n,matcoloring->ncolors);CHKERRQ(ierr); 2882128a86cSHong Zhang ierr = MatSetType(Bt_dense,MATSEQDENSE);CHKERRQ(ierr); 2892128a86cSHong Zhang ierr = MatSeqDenseSetPreallocation(Bt_dense,PETSC_NULL);CHKERRQ(ierr); 2902128a86cSHong Zhang Bt_dense->assembled = PETSC_TRUE; 2912128a86cSHong Zhang multtrans->Bt_den = Bt_dense; 2922128a86cSHong Zhang 2932128a86cSHong Zhang ierr = MatCreate(PETSC_COMM_SELF,&C_dense);CHKERRQ(ierr); 2942128a86cSHong Zhang ierr = MatSetSizes(C_dense,A->rmap->n,matcoloring->ncolors,A->rmap->n,matcoloring->ncolors);CHKERRQ(ierr); 2952128a86cSHong Zhang ierr = MatSetType(C_dense,MATSEQDENSE);CHKERRQ(ierr); 2962128a86cSHong Zhang ierr = MatSeqDenseSetPreallocation(C_dense,PETSC_NULL);CHKERRQ(ierr); 2972128a86cSHong Zhang Bt_dense->assembled = PETSC_TRUE; 2982128a86cSHong Zhang multtrans->ABt_den = C_dense; 299e8354b3eSHong Zhang ierr = PetscGetTime(&tf);CHKERRQ(ierr); 300e8354b3eSHong Zhang etime1 += tf - t0; 301f94ccd6cSHong Zhang 302f94ccd6cSHong Zhang #if defined(PETSC_USE_INFO) 303f94ccd6cSHong Zhang Mat_SeqAIJ *c=(Mat_SeqAIJ*)(*C)->data; 304f94ccd6cSHong 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)); 305f94ccd6cSHong Zhang ierr = PetscInfo5(*C,"Sym = GetColor %g + ColorCreate %g + MatDenseCreate %g + non-colorSym %g = %g\n",etime0,etime01,etime1,etime2,etime0+etime01+etime1+etime2); 306f94ccd6cSHong Zhang #endif 3072128a86cSHong Zhang } 308*e99be685SHong Zhang /* clean up */ 309*e99be685SHong Zhang ierr = MatDestroy(&Bt);CHKERRQ(ierr); 310*e99be685SHong Zhang ierr = MatRestoreSymbolicTranspose_SeqAIJ(B,&bti,&btj);CHKERRQ(ierr); 3112128a86cSHong Zhang 312f94ccd6cSHong Zhang 313f94ccd6cSHong Zhang 31481d82fe4SHong Zhang #if defined(INEFFICIENT_ALGORITHM) 31581d82fe4SHong Zhang /* The algorithm below computes am*bm sparse inner-product - inefficient! It will be deleted later. */ 3161fa1209cSHong Zhang PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 3171fa1209cSHong Zhang Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c; 3181fa1209cSHong Zhang PetscInt *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*ci,*cj,*acol,*bcol; 3191fa1209cSHong Zhang PetscInt am=A->rmap->N,bm=B->rmap->N; 3201fa1209cSHong Zhang PetscInt i,j,anzi,bnzj,cnzi,nlnk,*lnk,nspacedouble=0,ka,kb,index[1]; 3211fa1209cSHong Zhang MatScalar *ca; 3221fa1209cSHong Zhang PetscBT lnkbt; 3231fa1209cSHong Zhang PetscReal afill; 3245df89d91SHong Zhang 3251fa1209cSHong Zhang /* Allocate row pointer array ci */ 3261fa1209cSHong Zhang ierr = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr); 3271fa1209cSHong Zhang ci[0] = 0; 3281fa1209cSHong Zhang 3291fa1209cSHong Zhang /* Create and initialize a linked list for C columns */ 3301fa1209cSHong Zhang nlnk = bm+1; 3311fa1209cSHong Zhang ierr = PetscLLCreate(bm,bm,nlnk,lnk,lnkbt);CHKERRQ(ierr); 3321fa1209cSHong Zhang 3331fa1209cSHong Zhang /* Initial FreeSpace with size fill*(nnz(A)+nnz(B)) */ 3341fa1209cSHong Zhang ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr); 3351fa1209cSHong Zhang current_space = free_space; 3361fa1209cSHong Zhang 3371fa1209cSHong Zhang /* Determine symbolic info for each row of the product A*B^T: */ 3381fa1209cSHong Zhang for (i=0; i<am; i++) { 3391fa1209cSHong Zhang anzi = ai[i+1] - ai[i]; 3401fa1209cSHong Zhang cnzi = 0; 3411fa1209cSHong Zhang acol = aj + ai[i]; 3421fa1209cSHong Zhang for (j=0; j<bm; j++){ 3431fa1209cSHong Zhang bnzj = bi[j+1] - bi[j]; 3441fa1209cSHong Zhang bcol= bj + bi[j]; 34581d82fe4SHong Zhang /* sparse inner-product c(i,j)=A[i,:]*B[j,:]^T */ 3461fa1209cSHong Zhang ka = 0; kb = 0; 3471fa1209cSHong Zhang while (ka < anzi && kb < bnzj){ 34881d82fe4SHong Zhang while (acol[ka] < bcol[kb] && ka < anzi) ka++; 34981d82fe4SHong Zhang if (ka == anzi) break; 35081d82fe4SHong Zhang while (acol[ka] > bcol[kb] && kb < bnzj) kb++; 35181d82fe4SHong Zhang if (kb == bnzj) break; 3521fa1209cSHong Zhang if (acol[ka] == bcol[kb]){ /* add nonzero c(i,j) to lnk */ 3531fa1209cSHong Zhang index[0] = j; 3541fa1209cSHong Zhang ierr = PetscLLAdd(1,index,bm,nlnk,lnk,lnkbt);CHKERRQ(ierr); 3551fa1209cSHong Zhang cnzi++; 3561fa1209cSHong Zhang break; 3571fa1209cSHong Zhang } 3581fa1209cSHong Zhang } 3591fa1209cSHong Zhang } 3601fa1209cSHong Zhang 3611fa1209cSHong Zhang /* If free space is not available, make more free space */ 3621fa1209cSHong Zhang /* Double the amount of total space in the list */ 3631fa1209cSHong Zhang if (current_space->local_remaining<cnzi) { 3641fa1209cSHong Zhang ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 3651fa1209cSHong Zhang nspacedouble++; 3661fa1209cSHong Zhang } 3671fa1209cSHong Zhang 3681fa1209cSHong Zhang /* Copy data into free space, then initialize lnk */ 3691fa1209cSHong Zhang ierr = PetscLLClean(bm,bm,cnzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 3701fa1209cSHong Zhang current_space->array += cnzi; 3711fa1209cSHong Zhang current_space->local_used += cnzi; 3721fa1209cSHong Zhang current_space->local_remaining -= cnzi; 3731fa1209cSHong Zhang 3741fa1209cSHong Zhang ci[i+1] = ci[i] + cnzi; 3751fa1209cSHong Zhang } 3761fa1209cSHong Zhang 3771fa1209cSHong Zhang 3781fa1209cSHong Zhang /* Column indices are in the list of free space. 3791fa1209cSHong Zhang Allocate array cj, copy column indices to cj, and destroy list of free space */ 3801fa1209cSHong Zhang ierr = PetscMalloc((ci[am]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr); 3811fa1209cSHong Zhang ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr); 3821fa1209cSHong Zhang ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 3831fa1209cSHong Zhang 3841fa1209cSHong Zhang /* Allocate space for ca */ 3851fa1209cSHong Zhang ierr = PetscMalloc((ci[am]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr); 3861fa1209cSHong Zhang ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr); 3871fa1209cSHong Zhang 3881fa1209cSHong Zhang /* put together the new symbolic matrix */ 3891fa1209cSHong Zhang ierr = MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,am,bm,ci,cj,ca,C);CHKERRQ(ierr); 3901fa1209cSHong Zhang 3911fa1209cSHong Zhang /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 3921fa1209cSHong Zhang /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */ 3931fa1209cSHong Zhang c = (Mat_SeqAIJ *)((*C)->data); 3941fa1209cSHong Zhang c->free_a = PETSC_TRUE; 3951fa1209cSHong Zhang c->free_ij = PETSC_TRUE; 3961fa1209cSHong Zhang c->nonew = 0; 3971fa1209cSHong Zhang 3981fa1209cSHong Zhang /* set MatInfo */ 3991fa1209cSHong Zhang afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5; 4001fa1209cSHong Zhang if (afill < 1.0) afill = 1.0; 4011fa1209cSHong Zhang c->maxnz = ci[am]; 4021fa1209cSHong Zhang c->nz = ci[am]; 4031fa1209cSHong Zhang (*C)->info.mallocs = nspacedouble; 4041fa1209cSHong Zhang (*C)->info.fill_ratio_given = fill; 4051fa1209cSHong Zhang (*C)->info.fill_ratio_needed = afill; 4061fa1209cSHong Zhang 4071fa1209cSHong Zhang #if defined(PETSC_USE_INFO) 4081fa1209cSHong Zhang if (ci[am]) { 4091fa1209cSHong Zhang ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr); 4106fc122caSHong Zhang ierr = PetscInfo1((*C),"Use MatMatTransposeMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr); 4111fa1209cSHong Zhang } else { 4121fa1209cSHong Zhang ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr); 4131fa1209cSHong Zhang } 4141fa1209cSHong Zhang #endif 41581d82fe4SHong Zhang #endif 4165df89d91SHong Zhang PetscFunctionReturn(0); 4175df89d91SHong Zhang } 4185df89d91SHong Zhang 4196973769bSHong Zhang /* #define USE_ARRAY - for sparse dot product. Slower than !USE_ARRAY */ 4205df89d91SHong Zhang #undef __FUNCT__ 4216fc122caSHong Zhang #define __FUNCT__ "MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ" 4226fc122caSHong Zhang PetscErrorCode MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C) 4235df89d91SHong Zhang { 4245df89d91SHong Zhang PetscErrorCode ierr; 4255df89d91SHong Zhang Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data; 426e2cac8caSJed Brown PetscInt *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,anzi,bnzj,nexta,nextb,*acol,*bcol,brow; 42781d82fe4SHong Zhang PetscInt cm=C->rmap->n,*ci=c->i,*cj=c->j,i,j,cnzi,*ccol; 4285df89d91SHong Zhang PetscLogDouble flops=0.0; 4296973769bSHong Zhang MatScalar *aa=a->a,*aval,*ba=b->a,*bval,*ca=c->a,*cval; 4302128a86cSHong Zhang Mat_MatMatMultTrans *multtrans; 4312128a86cSHong Zhang PetscContainer container; 4326973769bSHong Zhang #if defined(USE_ARRAY) 4336973769bSHong Zhang MatScalar *spdot; 4346973769bSHong Zhang #endif 4355df89d91SHong Zhang 4365df89d91SHong Zhang PetscFunctionBegin; 4372128a86cSHong Zhang ierr = PetscObjectQuery((PetscObject)C,"Mat_MatMatMultTrans",(PetscObject *)&container);CHKERRQ(ierr); 4382128a86cSHong Zhang if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit"); 4392128a86cSHong Zhang ierr = PetscContainerGetPointer(container,(void **)&multtrans);CHKERRQ(ierr); 4402128a86cSHong Zhang if (multtrans->usecoloring){ 441b9af6bddSHong Zhang MatTransposeColoring matcoloring = multtrans->matcoloring; 442c8db22f6SHong Zhang Mat Bt_dense; 443c8db22f6SHong Zhang PetscInt m,n; 444b2d2b04fSHong Zhang PetscLogDouble t0,tf,etime0=0.0,etime1=0.0,etime2=0.0; 445a0b95ffbSSatish Balay Mat C_dense = multtrans->ABt_den; 446a0b95ffbSSatish Balay 4472128a86cSHong Zhang Bt_dense = multtrans->Bt_den; 448c8db22f6SHong Zhang ierr = MatGetLocalSize(Bt_dense,&m,&n);CHKERRQ(ierr); 449c8db22f6SHong Zhang 450b9af6bddSHong Zhang /* Get Bt_dense by Apply MatTransposeColoring to B */ 4512128a86cSHong Zhang ierr = PetscGetTime(&t0);CHKERRQ(ierr); 452b9af6bddSHong Zhang ierr = MatTransColoringApplySpToDen(matcoloring,B,Bt_dense);CHKERRQ(ierr); 4532128a86cSHong Zhang ierr = PetscGetTime(&tf);CHKERRQ(ierr); 454b2d2b04fSHong Zhang etime0 += tf - t0; 455c8db22f6SHong Zhang 456c8db22f6SHong Zhang /* C_dense = A*Bt_dense */ 4572128a86cSHong Zhang ierr = PetscGetTime(&t0);CHKERRQ(ierr); 4582128a86cSHong Zhang ierr = MatMatMultNumeric_SeqAIJ_SeqDense(A,Bt_dense,C_dense);CHKERRQ(ierr); 4592128a86cSHong Zhang ierr = PetscGetTime(&tf);CHKERRQ(ierr); 4602128a86cSHong Zhang etime2 += tf - t0; 461c8db22f6SHong Zhang 4622128a86cSHong Zhang /* Recover C from C_dense */ 4632128a86cSHong Zhang ierr = PetscGetTime(&t0);CHKERRQ(ierr); 464b9af6bddSHong Zhang ierr = MatTransColoringApplyDenToSp(matcoloring,C_dense,C);CHKERRQ(ierr); 4652128a86cSHong Zhang ierr = PetscGetTime(&tf);CHKERRQ(ierr); 4662128a86cSHong Zhang etime1 += tf - t0; 467f94ccd6cSHong Zhang #if defined(PETSC_USE_INFO) 468f94ccd6cSHong Zhang ierr = PetscInfo4(C,"Num = ColoringApply: %g %g + Mult_sp_dense %g = %g\n",etime0,etime1,etime2,etime0+etime1+etime2); 469f94ccd6cSHong Zhang #endif 470c8db22f6SHong Zhang PetscFunctionReturn(0); 471c8db22f6SHong Zhang } 472c8db22f6SHong Zhang 4736973769bSHong Zhang #if defined(USE_ARRAY) 4746973769bSHong Zhang /* allocate an array for implementing sparse inner-product */ 475e2cac8caSJed Brown ierr = PetscMalloc((A->cmap->n+1)*sizeof(MatScalar),&spdot);CHKERRQ(ierr); 476e2cac8caSJed Brown ierr = PetscMemzero(spdot,(A->cmap->n+1)*sizeof(MatScalar));CHKERRQ(ierr); 4776973769bSHong Zhang #endif 4786973769bSHong Zhang 47981d82fe4SHong Zhang /* clear old values in C */ 48081d82fe4SHong Zhang ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr); 4815df89d91SHong Zhang 4821fa1209cSHong Zhang for (i=0; i<cm; i++) { 48381d82fe4SHong Zhang anzi = ai[i+1] - ai[i]; 48481d82fe4SHong Zhang acol = aj + ai[i]; 4856973769bSHong Zhang aval = aa + ai[i]; 4861fa1209cSHong Zhang cnzi = ci[i+1] - ci[i]; 4871fa1209cSHong Zhang ccol = cj + ci[i]; 4886973769bSHong Zhang cval = ca + ci[i]; 4891fa1209cSHong Zhang for (j=0; j<cnzi; j++){ 49081d82fe4SHong Zhang brow = ccol[j]; 49181d82fe4SHong Zhang bnzj = bi[brow+1] - bi[brow]; 49281d82fe4SHong Zhang bcol = bj + bi[brow]; 4936973769bSHong Zhang bval = ba + bi[brow]; 4946973769bSHong Zhang 49581d82fe4SHong Zhang /* perform sparse inner-product c(i,j)=A[i,:]*B[j,:]^T */ 4966973769bSHong Zhang #if defined(USE_ARRAY) 4976973769bSHong Zhang /* put ba to spdot array */ 4986973769bSHong Zhang for (nextb=0; nextb<bnzj; nextb++) spdot[bcol[nextb]] = bval[nextb]; 4996973769bSHong Zhang /* c(i,j)=A[i,:]*B[j,:]^T */ 5006973769bSHong Zhang for (nexta=0; nexta<anzi; nexta++){ 5016973769bSHong Zhang cval[j] += spdot[acol[nexta]]*aval[nexta]; 5026973769bSHong Zhang } 5036973769bSHong Zhang /* zero spdot array */ 5046973769bSHong Zhang for (nextb=0; nextb<bnzj; nextb++) spdot[bcol[nextb]] = 0.0; 5056973769bSHong Zhang #else 50681d82fe4SHong Zhang nexta = 0; nextb = 0; 50781d82fe4SHong Zhang while (nexta<anzi && nextb<bnzj){ 50881d82fe4SHong Zhang while (acol[nexta] < bcol[nextb] && nexta < anzi) nexta++; 50981d82fe4SHong Zhang if (nexta == anzi) break; 51081d82fe4SHong Zhang while (acol[nexta] > bcol[nextb] && nextb < bnzj) nextb++; 51181d82fe4SHong Zhang if (nextb == bnzj) break; 51281d82fe4SHong Zhang if (acol[nexta] == bcol[nextb]){ 5136973769bSHong Zhang cval[j] += aval[nexta]*bval[nextb]; 51481d82fe4SHong Zhang nexta++; nextb++; 51581d82fe4SHong Zhang flops += 2; 5161fa1209cSHong Zhang } 5171fa1209cSHong Zhang } 5186973769bSHong Zhang #endif 51981d82fe4SHong Zhang } 52081d82fe4SHong Zhang } 52181d82fe4SHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 52281d82fe4SHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 52381d82fe4SHong Zhang ierr = PetscLogFlops(flops);CHKERRQ(ierr); 5246973769bSHong Zhang #if defined(USE_ARRAY) 5256973769bSHong Zhang ierr = PetscFree(spdot); 5266973769bSHong Zhang #endif 5275df89d91SHong Zhang PetscFunctionReturn(0); 5285df89d91SHong Zhang } 5295df89d91SHong Zhang 5305df89d91SHong Zhang #undef __FUNCT__ 53175648e8dSHong Zhang #define __FUNCT__ "MatTransposeMatMult_SeqAIJ_SeqAIJ" 53275648e8dSHong Zhang PetscErrorCode MatTransposeMatMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) { 5335df89d91SHong Zhang PetscErrorCode ierr; 5345df89d91SHong Zhang 5355df89d91SHong Zhang PetscFunctionBegin; 5365df89d91SHong Zhang if (scall == MAT_INITIAL_MATRIX){ 53775648e8dSHong Zhang ierr = MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr); 5385df89d91SHong Zhang } 53975648e8dSHong Zhang ierr = MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr); 5405df89d91SHong Zhang PetscFunctionReturn(0); 5415df89d91SHong Zhang } 5425df89d91SHong Zhang 5435df89d91SHong Zhang #undef __FUNCT__ 54475648e8dSHong Zhang #define __FUNCT__ "MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ" 54575648e8dSHong Zhang PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C) 5465df89d91SHong Zhang { 547bc011b1eSHong Zhang PetscErrorCode ierr; 548bc011b1eSHong Zhang Mat At; 54938baddfdSBarry Smith PetscInt *ati,*atj; 550bc011b1eSHong Zhang 551bc011b1eSHong Zhang PetscFunctionBegin; 552bc011b1eSHong Zhang /* create symbolic At */ 553bc011b1eSHong Zhang ierr = MatGetSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr); 554d0f46423SBarry Smith ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,A->cmap->n,A->rmap->n,ati,atj,PETSC_NULL,&At);CHKERRQ(ierr); 555bc011b1eSHong Zhang 556bc011b1eSHong Zhang /* get symbolic C=At*B */ 557bc011b1eSHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(At,B,fill,C);CHKERRQ(ierr); 558bc011b1eSHong Zhang 559bc011b1eSHong Zhang /* clean up */ 5606bf464f9SBarry Smith ierr = MatDestroy(&At);CHKERRQ(ierr); 561bc011b1eSHong Zhang ierr = MatRestoreSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr); 562bc011b1eSHong Zhang PetscFunctionReturn(0); 563bc011b1eSHong Zhang } 564bc011b1eSHong Zhang 565bc011b1eSHong Zhang #undef __FUNCT__ 56675648e8dSHong Zhang #define __FUNCT__ "MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ" 56775648e8dSHong Zhang PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C) 568bc011b1eSHong Zhang { 569bc011b1eSHong Zhang PetscErrorCode ierr; 5700fbc74f4SHong Zhang Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data; 571d0f46423SBarry Smith PetscInt am=A->rmap->n,anzi,*ai=a->i,*aj=a->j,*bi=b->i,*bj,bnzi,nextb; 572d0f46423SBarry Smith PetscInt cm=C->rmap->n,*ci=c->i,*cj=c->j,crow,*cjj,i,j,k; 573d13dce4bSSatish Balay PetscLogDouble flops=0.0; 5740fbc74f4SHong Zhang MatScalar *aa=a->a,*ba,*ca=c->a,*caj; 575bc011b1eSHong Zhang 576bc011b1eSHong Zhang PetscFunctionBegin; 577bc011b1eSHong Zhang /* clear old values in C */ 578bc011b1eSHong Zhang ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr); 579bc011b1eSHong Zhang 580bc011b1eSHong Zhang /* compute A^T*B using outer product (A^T)[:,i]*B[i,:] */ 581bc011b1eSHong Zhang for (i=0;i<am;i++) { 582bc011b1eSHong Zhang bj = b->j + bi[i]; 583bc011b1eSHong Zhang ba = b->a + bi[i]; 584bc011b1eSHong Zhang bnzi = bi[i+1] - bi[i]; 585bc011b1eSHong Zhang anzi = ai[i+1] - ai[i]; 586bc011b1eSHong Zhang for (j=0; j<anzi; j++) { 587bc011b1eSHong Zhang nextb = 0; 5880fbc74f4SHong Zhang crow = *aj++; 589bc011b1eSHong Zhang cjj = cj + ci[crow]; 590bc011b1eSHong Zhang caj = ca + ci[crow]; 591bc011b1eSHong Zhang /* perform sparse axpy operation. Note cjj includes bj. */ 592bc011b1eSHong Zhang for (k=0; nextb<bnzi; k++) { 5930fbc74f4SHong Zhang if (cjj[k] == *(bj+nextb)) { /* ccol == bcol */ 5940fbc74f4SHong Zhang caj[k] += (*aa)*(*(ba+nextb)); 595bc011b1eSHong Zhang nextb++; 596bc011b1eSHong Zhang } 597bc011b1eSHong Zhang } 598bc011b1eSHong Zhang flops += 2*bnzi; 5990fbc74f4SHong Zhang aa++; 600bc011b1eSHong Zhang } 601bc011b1eSHong Zhang } 602bc011b1eSHong Zhang 603bc011b1eSHong Zhang /* Assemble the final matrix and clean up */ 604bc011b1eSHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 605bc011b1eSHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 606bc011b1eSHong Zhang ierr = PetscLogFlops(flops);CHKERRQ(ierr); 607bc011b1eSHong Zhang PetscFunctionReturn(0); 608bc011b1eSHong Zhang } 6099123193aSHong Zhang 610ec01deb9SMatthew Knepley EXTERN_C_BEGIN 6119123193aSHong Zhang #undef __FUNCT__ 6129123193aSHong Zhang #define __FUNCT__ "MatMatMult_SeqAIJ_SeqDense" 6139123193aSHong Zhang PetscErrorCode MatMatMult_SeqAIJ_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 6149123193aSHong Zhang { 6159123193aSHong Zhang PetscErrorCode ierr; 6169123193aSHong Zhang 6179123193aSHong Zhang PetscFunctionBegin; 6189123193aSHong Zhang if (scall == MAT_INITIAL_MATRIX){ 6199123193aSHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqDense(A,B,fill,C);CHKERRQ(ierr); 6209123193aSHong Zhang } 6219123193aSHong Zhang ierr = MatMatMultNumeric_SeqAIJ_SeqDense(A,B,*C);CHKERRQ(ierr); 6229123193aSHong Zhang PetscFunctionReturn(0); 6239123193aSHong Zhang } 624ec01deb9SMatthew Knepley EXTERN_C_END 625edd81eecSMatthew Knepley 6269123193aSHong Zhang #undef __FUNCT__ 6279123193aSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqDense" 6289123193aSHong Zhang PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C) 6299123193aSHong Zhang { 6309123193aSHong Zhang PetscErrorCode ierr; 6319123193aSHong Zhang 6329123193aSHong Zhang PetscFunctionBegin; 6335a586d82SBarry Smith ierr = MatMatMultSymbolic_SeqDense_SeqDense(A,B,0.0,C);CHKERRQ(ierr); 6349123193aSHong Zhang PetscFunctionReturn(0); 6359123193aSHong Zhang } 6369123193aSHong Zhang 6379123193aSHong Zhang #undef __FUNCT__ 6389123193aSHong Zhang #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqDense" 6399123193aSHong Zhang PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqDense(Mat A,Mat B,Mat C) 6409123193aSHong Zhang { 641f32d5d43SBarry Smith Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data; 6429123193aSHong Zhang PetscErrorCode ierr; 643dd6ea824SBarry Smith PetscScalar *b,*c,r1,r2,r3,r4,*b1,*b2,*b3,*b4; 644dd6ea824SBarry Smith MatScalar *aa; 645d0f46423SBarry Smith PetscInt cm=C->rmap->n, cn=B->cmap->n, bm=B->rmap->n, col, i,j,n,*aj, am = A->rmap->n; 646f32d5d43SBarry Smith PetscInt am2 = 2*am, am3 = 3*am, bm4 = 4*bm,colam; 6479123193aSHong Zhang 6489123193aSHong Zhang PetscFunctionBegin; 649f32d5d43SBarry Smith if (!cm || !cn) PetscFunctionReturn(0); 650e32f2f54SBarry 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); 651e32f2f54SBarry 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); 652e32f2f54SBarry 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); 653f32d5d43SBarry Smith ierr = MatGetArray(B,&b);CHKERRQ(ierr); 654f32d5d43SBarry Smith ierr = MatGetArray(C,&c);CHKERRQ(ierr); 655f32d5d43SBarry Smith b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm; 656f32d5d43SBarry Smith for (col=0; col<cn-4; col += 4){ /* over columns of C */ 657f32d5d43SBarry Smith colam = col*am; 658f32d5d43SBarry Smith for (i=0; i<am; i++) { /* over rows of C in those columns */ 659f32d5d43SBarry Smith r1 = r2 = r3 = r4 = 0.0; 660f32d5d43SBarry Smith n = a->i[i+1] - a->i[i]; 661f32d5d43SBarry Smith aj = a->j + a->i[i]; 662f32d5d43SBarry Smith aa = a->a + a->i[i]; 663f32d5d43SBarry Smith for (j=0; j<n; j++) { 664f32d5d43SBarry Smith r1 += (*aa)*b1[*aj]; 665f32d5d43SBarry Smith r2 += (*aa)*b2[*aj]; 666f32d5d43SBarry Smith r3 += (*aa)*b3[*aj]; 667f32d5d43SBarry Smith r4 += (*aa++)*b4[*aj++]; 6689123193aSHong Zhang } 669f32d5d43SBarry Smith c[colam + i] = r1; 670f32d5d43SBarry Smith c[colam + am + i] = r2; 671f32d5d43SBarry Smith c[colam + am2 + i] = r3; 672f32d5d43SBarry Smith c[colam + am3 + i] = r4; 673f32d5d43SBarry Smith } 674f32d5d43SBarry Smith b1 += bm4; 675f32d5d43SBarry Smith b2 += bm4; 676f32d5d43SBarry Smith b3 += bm4; 677f32d5d43SBarry Smith b4 += bm4; 678f32d5d43SBarry Smith } 679f32d5d43SBarry Smith for (;col<cn; col++){ /* over extra columns of C */ 680f32d5d43SBarry Smith for (i=0; i<am; i++) { /* over rows of C in those columns */ 681f32d5d43SBarry Smith r1 = 0.0; 682f32d5d43SBarry Smith n = a->i[i+1] - a->i[i]; 683f32d5d43SBarry Smith aj = a->j + a->i[i]; 684f32d5d43SBarry Smith aa = a->a + a->i[i]; 685f32d5d43SBarry Smith 686f32d5d43SBarry Smith for (j=0; j<n; j++) { 687f32d5d43SBarry Smith r1 += (*aa++)*b1[*aj++]; 688f32d5d43SBarry Smith } 689f32d5d43SBarry Smith c[col*am + i] = r1; 690f32d5d43SBarry Smith } 691f32d5d43SBarry Smith b1 += bm; 692f32d5d43SBarry Smith } 693dc0b31edSSatish Balay ierr = PetscLogFlops(cn*(2.0*a->nz));CHKERRQ(ierr); 694f32d5d43SBarry Smith ierr = MatRestoreArray(B,&b);CHKERRQ(ierr); 695f32d5d43SBarry Smith ierr = MatRestoreArray(C,&c);CHKERRQ(ierr); 696f32d5d43SBarry Smith ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 697f32d5d43SBarry Smith ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 698f32d5d43SBarry Smith PetscFunctionReturn(0); 699f32d5d43SBarry Smith } 700f32d5d43SBarry Smith 701f32d5d43SBarry Smith /* 7024324174eSBarry Smith Note very similar to MatMult_SeqAIJ(), should generate both codes from same base 703f32d5d43SBarry Smith */ 704f32d5d43SBarry Smith #undef __FUNCT__ 705f32d5d43SBarry Smith #define __FUNCT__ "MatMatMultNumericAdd_SeqAIJ_SeqDense" 706f32d5d43SBarry Smith PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat A,Mat B,Mat C) 707f32d5d43SBarry Smith { 708f32d5d43SBarry Smith Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data; 709f32d5d43SBarry Smith PetscErrorCode ierr; 710dd6ea824SBarry Smith PetscScalar *b,*c,r1,r2,r3,r4,*b1,*b2,*b3,*b4; 711dd6ea824SBarry Smith MatScalar *aa; 712d0f46423SBarry 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; 7134324174eSBarry Smith PetscInt am2 = 2*am, am3 = 3*am, bm4 = 4*bm,colam,*ridx; 714f32d5d43SBarry Smith 715f32d5d43SBarry Smith PetscFunctionBegin; 716f32d5d43SBarry Smith if (!cm || !cn) PetscFunctionReturn(0); 717f32d5d43SBarry Smith ierr = MatGetArray(B,&b);CHKERRQ(ierr); 718f32d5d43SBarry Smith ierr = MatGetArray(C,&c);CHKERRQ(ierr); 719f32d5d43SBarry Smith b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm; 7204324174eSBarry Smith 7214324174eSBarry Smith if (a->compressedrow.use){ /* use compressed row format */ 7224324174eSBarry Smith for (col=0; col<cn-4; col += 4){ /* over columns of C */ 7234324174eSBarry Smith colam = col*am; 7244324174eSBarry Smith arm = a->compressedrow.nrows; 7254324174eSBarry Smith ii = a->compressedrow.i; 7264324174eSBarry Smith ridx = a->compressedrow.rindex; 7274324174eSBarry Smith for (i=0; i<arm; i++) { /* over rows of C in those columns */ 7284324174eSBarry Smith r1 = r2 = r3 = r4 = 0.0; 7294324174eSBarry Smith n = ii[i+1] - ii[i]; 7304324174eSBarry Smith aj = a->j + ii[i]; 7314324174eSBarry Smith aa = a->a + ii[i]; 7324324174eSBarry Smith for (j=0; j<n; j++) { 7334324174eSBarry Smith r1 += (*aa)*b1[*aj]; 7344324174eSBarry Smith r2 += (*aa)*b2[*aj]; 7354324174eSBarry Smith r3 += (*aa)*b3[*aj]; 7364324174eSBarry Smith r4 += (*aa++)*b4[*aj++]; 7374324174eSBarry Smith } 7384324174eSBarry Smith c[colam + ridx[i]] += r1; 7394324174eSBarry Smith c[colam + am + ridx[i]] += r2; 7404324174eSBarry Smith c[colam + am2 + ridx[i]] += r3; 7414324174eSBarry Smith c[colam + am3 + ridx[i]] += r4; 7424324174eSBarry Smith } 7434324174eSBarry Smith b1 += bm4; 7444324174eSBarry Smith b2 += bm4; 7454324174eSBarry Smith b3 += bm4; 7464324174eSBarry Smith b4 += bm4; 7474324174eSBarry Smith } 7484324174eSBarry Smith for (;col<cn; col++){ /* over extra columns of C */ 7494324174eSBarry Smith colam = col*am; 7504324174eSBarry Smith arm = a->compressedrow.nrows; 7514324174eSBarry Smith ii = a->compressedrow.i; 7524324174eSBarry Smith ridx = a->compressedrow.rindex; 7534324174eSBarry Smith for (i=0; i<arm; i++) { /* over rows of C in those columns */ 7544324174eSBarry Smith r1 = 0.0; 7554324174eSBarry Smith n = ii[i+1] - ii[i]; 7564324174eSBarry Smith aj = a->j + ii[i]; 7574324174eSBarry Smith aa = a->a + ii[i]; 7584324174eSBarry Smith 7594324174eSBarry Smith for (j=0; j<n; j++) { 7604324174eSBarry Smith r1 += (*aa++)*b1[*aj++]; 7614324174eSBarry Smith } 7624324174eSBarry Smith c[col*am + ridx[i]] += r1; 7634324174eSBarry Smith } 7644324174eSBarry Smith b1 += bm; 7654324174eSBarry Smith } 7664324174eSBarry Smith } else { 767f32d5d43SBarry Smith for (col=0; col<cn-4; col += 4){ /* over columns of C */ 768f32d5d43SBarry Smith colam = col*am; 769f32d5d43SBarry Smith for (i=0; i<am; i++) { /* over rows of C in those columns */ 770f32d5d43SBarry Smith r1 = r2 = r3 = r4 = 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 for (j=0; j<n; j++) { 775f32d5d43SBarry Smith r1 += (*aa)*b1[*aj]; 776f32d5d43SBarry Smith r2 += (*aa)*b2[*aj]; 777f32d5d43SBarry Smith r3 += (*aa)*b3[*aj]; 778f32d5d43SBarry Smith r4 += (*aa++)*b4[*aj++]; 779f32d5d43SBarry Smith } 780f32d5d43SBarry Smith c[colam + i] += r1; 781f32d5d43SBarry Smith c[colam + am + i] += r2; 782f32d5d43SBarry Smith c[colam + am2 + i] += r3; 783f32d5d43SBarry Smith c[colam + am3 + i] += r4; 784f32d5d43SBarry Smith } 785f32d5d43SBarry Smith b1 += bm4; 786f32d5d43SBarry Smith b2 += bm4; 787f32d5d43SBarry Smith b3 += bm4; 788f32d5d43SBarry Smith b4 += bm4; 789f32d5d43SBarry Smith } 790f32d5d43SBarry Smith for (;col<cn; col++){ /* over extra columns of C */ 791f32d5d43SBarry Smith for (i=0; i<am; i++) { /* over rows of C in those columns */ 792f32d5d43SBarry Smith r1 = 0.0; 793f32d5d43SBarry Smith n = a->i[i+1] - a->i[i]; 794f32d5d43SBarry Smith aj = a->j + a->i[i]; 795f32d5d43SBarry Smith aa = a->a + a->i[i]; 796f32d5d43SBarry Smith 797f32d5d43SBarry Smith for (j=0; j<n; j++) { 798f32d5d43SBarry Smith r1 += (*aa++)*b1[*aj++]; 799f32d5d43SBarry Smith } 800f32d5d43SBarry Smith c[col*am + i] += r1; 801f32d5d43SBarry Smith } 802f32d5d43SBarry Smith b1 += bm; 803f32d5d43SBarry Smith } 8044324174eSBarry Smith } 805dc0b31edSSatish Balay ierr = PetscLogFlops(cn*2.0*a->nz);CHKERRQ(ierr); 806f32d5d43SBarry Smith ierr = MatRestoreArray(B,&b);CHKERRQ(ierr); 807f32d5d43SBarry Smith ierr = MatRestoreArray(C,&c);CHKERRQ(ierr); 8089123193aSHong Zhang PetscFunctionReturn(0); 8099123193aSHong Zhang } 810b1683b59SHong Zhang 811b1683b59SHong Zhang #undef __FUNCT__ 812b9af6bddSHong Zhang #define __FUNCT__ "MatTransColoringApplySpToDen_SeqAIJ" 813b9af6bddSHong Zhang PetscErrorCode MatTransColoringApplySpToDen_SeqAIJ(MatTransposeColoring coloring,Mat B,Mat Btdense) 814c8db22f6SHong Zhang { 815c8db22f6SHong Zhang PetscErrorCode ierr; 8162f5f1f90SHong Zhang Mat_SeqAIJ *b = (Mat_SeqAIJ*)B->data; 8172f5f1f90SHong Zhang Mat_SeqDense *btdense = (Mat_SeqDense*)Btdense->data; 8182f5f1f90SHong Zhang PetscInt *bi=b->i,*bj=b->j; 8192f5f1f90SHong Zhang PetscInt m=Btdense->rmap->n,n=Btdense->cmap->n,j,k,l,col,anz,*btcol,brow,ncolumns; 8202f5f1f90SHong Zhang MatScalar *btval,*btval_den,*ba=b->a; 8212f5f1f90SHong Zhang PetscInt *columns=coloring->columns,*colorforcol=coloring->colorforcol,ncolors=coloring->ncolors; 822c8db22f6SHong Zhang 823c8db22f6SHong Zhang PetscFunctionBegin; 8242f5f1f90SHong Zhang btval_den=btdense->v; 8252f5f1f90SHong Zhang ierr = PetscMemzero(btval_den,(m*n)*sizeof(MatScalar));CHKERRQ(ierr); 8262f5f1f90SHong Zhang for (k=0; k<ncolors; k++) { 8272f5f1f90SHong Zhang ncolumns = coloring->ncolumns[k]; 8282f5f1f90SHong Zhang for (l=0; l<ncolumns; l++) { /* insert a row of B to a column of Btdense */ 829d2853540SHong Zhang col = *(columns + colorforcol[k] + l); 8302f5f1f90SHong Zhang btcol = bj + bi[col]; 8312f5f1f90SHong Zhang btval = ba + bi[col]; 8322f5f1f90SHong Zhang anz = bi[col+1] - bi[col]; 833d2853540SHong Zhang for (j=0; j<anz; j++){ 8342f5f1f90SHong Zhang brow = btcol[j]; 8352f5f1f90SHong Zhang btval_den[brow] = btval[j]; 836c8db22f6SHong Zhang } 837c8db22f6SHong Zhang } 8382f5f1f90SHong Zhang btval_den += m; 839c8db22f6SHong Zhang } 840c8db22f6SHong Zhang PetscFunctionReturn(0); 841c8db22f6SHong Zhang } 842c8db22f6SHong Zhang 843c8db22f6SHong Zhang #undef __FUNCT__ 844b9af6bddSHong Zhang #define __FUNCT__ "MatTransColoringApplyDenToSp_SeqAIJ" 845b9af6bddSHong Zhang PetscErrorCode MatTransColoringApplyDenToSp_SeqAIJ(MatTransposeColoring matcoloring,Mat Cden,Mat Csp) 8468972f759SHong Zhang { 8478972f759SHong Zhang PetscErrorCode ierr; 848b2d2b04fSHong Zhang Mat_SeqAIJ *csp = (Mat_SeqAIJ*)Csp->data; 8492f5f1f90SHong Zhang PetscInt k,l,*row,*idx,m,ncolors=matcoloring->ncolors,nrows; 850b2d2b04fSHong Zhang PetscScalar *ca_den,*cp_den,*ca=csp->a; 8512f5f1f90SHong Zhang PetscInt *rows=matcoloring->rows,*spidx=matcoloring->columnsforspidx,*colorforrow=matcoloring->colorforrow; 8528972f759SHong Zhang 8538972f759SHong Zhang PetscFunctionBegin; 8548972f759SHong Zhang ierr = MatGetLocalSize(Csp,&m,PETSC_NULL);CHKERRQ(ierr); 855b2d2b04fSHong Zhang ierr = MatGetArray(Cden,&ca_den);CHKERRQ(ierr); 856b2d2b04fSHong Zhang cp_den = ca_den; 8572f5f1f90SHong Zhang for (k=0; k<ncolors; k++) { 8582f5f1f90SHong Zhang nrows = matcoloring->nrows[k]; 8592f5f1f90SHong Zhang row = rows + colorforrow[k]; 8602f5f1f90SHong Zhang idx = spidx + colorforrow[k]; 8612f5f1f90SHong Zhang for (l=0; l<nrows; l++){ 8622f5f1f90SHong Zhang ca[idx[l]] = cp_den[row[l]]; 863b2d2b04fSHong Zhang } 864b2d2b04fSHong Zhang cp_den += m; 865b2d2b04fSHong Zhang } 866b2d2b04fSHong Zhang ierr = MatRestoreArray(Cden,&ca_den);CHKERRQ(ierr); 8678972f759SHong Zhang PetscFunctionReturn(0); 8688972f759SHong Zhang } 8698972f759SHong Zhang 870*e99be685SHong Zhang /* 871*e99be685SHong Zhang MatGetColumnIJ_SeqAIJ_Color() and MatRestoreColumnIJ_SeqAIJ_Color() are customized from 872*e99be685SHong Zhang MatGetColumnIJ_SeqAIJ() and MatRestoreColumnIJ_SeqAIJ() by adding an output 873*e99be685SHong Zhang spidx[], index of a->j, to be used for setting 'columnsforspidx' in MatTransposeColoringCreate_SeqAIJ(). 874*e99be685SHong Zhang */ 875*e99be685SHong Zhang #undef __FUNCT__ 876*e99be685SHong Zhang #define __FUNCT__ "MatGetColumnIJ_SeqAIJ_Color" 877*e99be685SHong Zhang PetscErrorCode MatGetColumnIJ_SeqAIJ_Color(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool inodecompressed,PetscInt *nn,PetscInt *ia[],PetscInt *ja[],PetscInt *spidx[],PetscBool *done) 878*e99be685SHong Zhang { 879*e99be685SHong Zhang Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 880*e99be685SHong Zhang PetscErrorCode ierr; 881*e99be685SHong Zhang PetscInt i,*collengths,*cia,*cja,n = A->cmap->n,m = A->rmap->n; 882*e99be685SHong Zhang PetscInt nz = a->i[m],row,*jj,mr,col; 883*e99be685SHong Zhang PetscInt *cspidx; 884*e99be685SHong Zhang 885*e99be685SHong Zhang PetscFunctionBegin; 886*e99be685SHong Zhang *nn = n; 887*e99be685SHong Zhang if (!ia) PetscFunctionReturn(0); 888*e99be685SHong Zhang if (symmetric) { 889*e99be685SHong Zhang SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"MatGetColumnIJ_SeqAIJ_Color() not supported for the case symmetric"); 890*e99be685SHong Zhang ierr = MatToSymmetricIJ_SeqAIJ(A->rmap->n,a->i,a->j,0,oshift,ia,ja);CHKERRQ(ierr); 891*e99be685SHong Zhang } else { 892*e99be685SHong Zhang ierr = PetscMalloc((n+1)*sizeof(PetscInt),&collengths);CHKERRQ(ierr); 893*e99be685SHong Zhang ierr = PetscMemzero(collengths,n*sizeof(PetscInt));CHKERRQ(ierr); 894*e99be685SHong Zhang ierr = PetscMalloc((n+1)*sizeof(PetscInt),&cia);CHKERRQ(ierr); 895*e99be685SHong Zhang ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&cja);CHKERRQ(ierr); 896*e99be685SHong Zhang ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&cspidx);CHKERRQ(ierr); 897*e99be685SHong Zhang jj = a->j; 898*e99be685SHong Zhang for (i=0; i<nz; i++) { 899*e99be685SHong Zhang collengths[jj[i]]++; 900*e99be685SHong Zhang } 901*e99be685SHong Zhang cia[0] = oshift; 902*e99be685SHong Zhang for (i=0; i<n; i++) { 903*e99be685SHong Zhang cia[i+1] = cia[i] + collengths[i]; 904*e99be685SHong Zhang } 905*e99be685SHong Zhang ierr = PetscMemzero(collengths,n*sizeof(PetscInt));CHKERRQ(ierr); 906*e99be685SHong Zhang jj = a->j; 907*e99be685SHong Zhang for (row=0; row<m; row++) { 908*e99be685SHong Zhang mr = a->i[row+1] - a->i[row]; 909*e99be685SHong Zhang for (i=0; i<mr; i++) { 910*e99be685SHong Zhang col = *jj++; 911*e99be685SHong Zhang cspidx[cia[col] + collengths[col] - oshift] = a->i[row] + i; /* index of a->j */ 912*e99be685SHong Zhang cja[cia[col] + collengths[col]++ - oshift] = row + oshift; 913*e99be685SHong Zhang } 914*e99be685SHong Zhang } 915*e99be685SHong Zhang ierr = PetscFree(collengths);CHKERRQ(ierr); 916*e99be685SHong Zhang *ia = cia; *ja = cja; 917*e99be685SHong Zhang *spidx = cspidx; 918*e99be685SHong Zhang } 919*e99be685SHong Zhang PetscFunctionReturn(0); 920*e99be685SHong Zhang } 921*e99be685SHong Zhang 922*e99be685SHong Zhang #undef __FUNCT__ 923*e99be685SHong Zhang #define __FUNCT__ "MatRestoreColumnIJ_SeqAIJ_Color" 924*e99be685SHong Zhang PetscErrorCode MatRestoreColumnIJ_SeqAIJ_Color(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool inodecompressed,PetscInt *n,PetscInt *ia[],PetscInt *ja[],PetscInt *spidx[],PetscBool *done) 925*e99be685SHong Zhang { 926*e99be685SHong Zhang PetscErrorCode ierr; 927*e99be685SHong Zhang 928*e99be685SHong Zhang PetscFunctionBegin; 929*e99be685SHong Zhang if (!ia) PetscFunctionReturn(0); 930*e99be685SHong Zhang 931*e99be685SHong Zhang ierr = PetscFree(*ia);CHKERRQ(ierr); 932*e99be685SHong Zhang ierr = PetscFree(*ja);CHKERRQ(ierr); 933*e99be685SHong Zhang ierr = PetscFree(*spidx);CHKERRQ(ierr); 934*e99be685SHong Zhang PetscFunctionReturn(0); 935*e99be685SHong Zhang } 936*e99be685SHong Zhang 9378972f759SHong Zhang #undef __FUNCT__ 938b9af6bddSHong Zhang #define __FUNCT__ "MatTransposeColoringCreate_SeqAIJ" 939b9af6bddSHong Zhang PetscErrorCode MatTransposeColoringCreate_SeqAIJ(Mat mat,ISColoring iscoloring,MatTransposeColoring c) 940b1683b59SHong Zhang { 941b1683b59SHong Zhang PetscErrorCode ierr; 942d2853540SHong Zhang PetscInt i,n,nrows,N,j,k,m,*row_idx,*ci,*cj,ncols,col,cm; 943b1683b59SHong Zhang const PetscInt *is; 944b28d80bdSHong Zhang PetscInt nis = iscoloring->n,*rowhit,bs = 1; 945b1683b59SHong Zhang IS *isa; 946b1683b59SHong Zhang PetscBool done; 947b1683b59SHong Zhang PetscBool flg1,flg2; 948b28d80bdSHong Zhang Mat_SeqAIJ *csp = (Mat_SeqAIJ*)mat->data; 949*e99be685SHong Zhang PetscInt *colorforrow,*rows,*rows_i,*columnsforspidx,*columnsforspidx_i,*idxhit,*spidx; 950d2853540SHong Zhang PetscInt *colorforcol,*columns,*columns_i; 951b1683b59SHong Zhang 952b1683b59SHong Zhang PetscFunctionBegin; 953b1683b59SHong Zhang ierr = ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);CHKERRQ(ierr); 954*e99be685SHong Zhang 955b1683b59SHong Zhang /* this is ugly way to get blocksize but cannot call MatGetBlockSize() because AIJ can have bs > 1 */ 956b1683b59SHong Zhang ierr = PetscTypeCompare((PetscObject)mat,MATSEQBAIJ,&flg1);CHKERRQ(ierr); 957b1683b59SHong Zhang ierr = PetscTypeCompare((PetscObject)mat,MATMPIBAIJ,&flg2);CHKERRQ(ierr); 958b1683b59SHong Zhang if (flg1 || flg2) { 959b1683b59SHong Zhang ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr); 960b1683b59SHong Zhang } 961b1683b59SHong Zhang 962b1683b59SHong Zhang N = mat->cmap->N/bs; 963b1683b59SHong Zhang c->M = mat->rmap->N/bs; /* set total rows, columns and local rows */ 964b1683b59SHong Zhang c->N = mat->cmap->N/bs; 965b1683b59SHong Zhang c->m = mat->rmap->N/bs; 966b1683b59SHong Zhang c->rstart = 0; 967b1683b59SHong Zhang 968b1683b59SHong Zhang c->ncolors = nis; 969b1683b59SHong Zhang ierr = PetscMalloc(nis*sizeof(PetscInt),&c->ncolumns);CHKERRQ(ierr); 970b28d80bdSHong Zhang ierr = PetscMalloc(nis*sizeof(PetscInt),&c->nrows);CHKERRQ(ierr); 971d2853540SHong Zhang ierr = PetscMalloc2(csp->nz+1,PetscInt,&rows,csp->nz+1,PetscInt,&columnsforspidx);CHKERRQ(ierr); 972d2853540SHong Zhang ierr = PetscMalloc((nis+1)*sizeof(PetscInt),&colorforrow);CHKERRQ(ierr); 973d2853540SHong Zhang colorforrow[0] = 0; 974d2853540SHong Zhang rows_i = rows; 975d2853540SHong Zhang columnsforspidx_i = columnsforspidx; 976b1683b59SHong Zhang 977d2853540SHong Zhang ierr = PetscMalloc((nis+1)*sizeof(PetscInt),&colorforcol);CHKERRQ(ierr); 978d2853540SHong Zhang ierr = PetscMalloc((N+1)*sizeof(PetscInt),&columns);CHKERRQ(ierr); 979d2853540SHong Zhang colorforcol[0] = 0; 980d2853540SHong Zhang columns_i = columns; 981d2853540SHong Zhang 982*e99be685SHong Zhang ierr = MatGetColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,&done);CHKERRQ(ierr); /* column-wise storage of mat */ 983b1683b59SHong Zhang if (!done) SETERRQ1(((PetscObject)mat)->comm,PETSC_ERR_SUP,"MatGetColumnIJ() not supported for matrix type %s",((PetscObject)mat)->type_name); 984b1683b59SHong Zhang 985b28d80bdSHong Zhang cm = c->m; 986b28d80bdSHong Zhang ierr = PetscMalloc((cm+1)*sizeof(PetscInt),&rowhit);CHKERRQ(ierr); 987*e99be685SHong Zhang ierr = PetscMalloc((cm+1)*sizeof(PetscInt),&idxhit);CHKERRQ(ierr); 988b1683b59SHong Zhang for (i=0; i<nis; i++) { 989b1683b59SHong Zhang ierr = ISGetLocalSize(isa[i],&n);CHKERRQ(ierr); 990b1683b59SHong Zhang ierr = ISGetIndices(isa[i],&is);CHKERRQ(ierr); 991b1683b59SHong Zhang c->ncolumns[i] = n; 992b1683b59SHong Zhang if (n) { 993d2853540SHong Zhang ierr = PetscMemcpy(columns_i,is,n*sizeof(PetscInt));CHKERRQ(ierr); 994b1683b59SHong Zhang } 995d2853540SHong Zhang colorforcol[i+1] = colorforcol[i] + n; 996d2853540SHong Zhang columns_i += n; 997b1683b59SHong Zhang 998b1683b59SHong Zhang /* fast, crude version requires O(N*N) work */ 999e8354b3eSHong Zhang ierr = PetscMemzero(rowhit,cm*sizeof(PetscInt));CHKERRQ(ierr); 1000*e99be685SHong Zhang 1001b1683b59SHong Zhang /* loop over columns*/ 1002b1683b59SHong Zhang for (j=0; j<n; j++) { 1003b1683b59SHong Zhang col = is[j]; 1004d2853540SHong Zhang row_idx = cj + ci[col]; 1005b1683b59SHong Zhang m = ci[col+1] - ci[col]; 1006b1683b59SHong Zhang /* loop over columns marking them in rowhit */ 1007b1683b59SHong Zhang for (k=0; k<m; k++) { 1008*e99be685SHong Zhang idxhit[*row_idx] = spidx[ci[col] + k]; 1009d2853540SHong Zhang rowhit[*row_idx++] = col + 1; 1010b1683b59SHong Zhang } 1011b1683b59SHong Zhang } 1012b1683b59SHong Zhang /* count the number of hits */ 1013b1683b59SHong Zhang nrows = 0; 1014e8354b3eSHong Zhang for (j=0; j<cm; j++) { 1015b1683b59SHong Zhang if (rowhit[j]) nrows++; 1016b1683b59SHong Zhang } 1017b1683b59SHong Zhang c->nrows[i] = nrows; 1018d2853540SHong Zhang colorforrow[i+1] = colorforrow[i] + nrows; 1019d2853540SHong Zhang 1020b1683b59SHong Zhang nrows = 0; 1021e8354b3eSHong Zhang for (j=0; j<cm; j++) { 1022b1683b59SHong Zhang if (rowhit[j]) { 1023d2853540SHong Zhang rows_i[nrows] = j; 1024*e99be685SHong Zhang columnsforspidx_i[nrows] = idxhit[j]; 1025b1683b59SHong Zhang nrows++; 1026b1683b59SHong Zhang } 1027b1683b59SHong Zhang } 1028b1683b59SHong Zhang ierr = ISRestoreIndices(isa[i],&is);CHKERRQ(ierr); 1029d2853540SHong Zhang rows_i += nrows; columnsforspidx_i += nrows; 1030b1683b59SHong Zhang } 1031*e99be685SHong Zhang ierr = MatRestoreColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,&done);CHKERRQ(ierr); 1032b28d80bdSHong Zhang ierr = PetscFree(rowhit);CHKERRQ(ierr); 1033b1683b59SHong Zhang ierr = ISColoringRestoreIS(iscoloring,&isa);CHKERRQ(ierr); 1034d2853540SHong Zhang #if defined(PETSC_USE_DEBUG) 1035d2853540SHong Zhang if (csp->nz != colorforrow[nis]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"csp->nz %d != colorforrow[nis] %d",csp->nz,colorforrow[nis]); 1036d2853540SHong Zhang #endif 1037b28d80bdSHong Zhang 1038d2853540SHong Zhang c->colorforrow = colorforrow; 1039d2853540SHong Zhang c->rows = rows; 1040d2853540SHong Zhang c->columnsforspidx = columnsforspidx; 1041d2853540SHong Zhang c->colorforcol = colorforcol; 1042d2853540SHong Zhang c->columns = columns; 1043f94ccd6cSHong Zhang ierr = PetscFree(idxhit);CHKERRQ(ierr); 1044b1683b59SHong Zhang PetscFunctionReturn(0); 1045b1683b59SHong Zhang } 1046