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*/ 11e9e4536cSHong Zhang /* 12e9e4536cSHong Zhang #define DEBUG_MATMATMULT 13e9e4536cSHong Zhang */ 14ec01deb9SMatthew Knepley EXTERN_C_BEGIN 156284ec50SHong Zhang #undef __FUNCT__ 165c66b693SKris Buschelman #define __FUNCT__ "MatMatMult_SeqAIJ_SeqAIJ" 1738baddfdSBarry Smith PetscErrorCode MatMatMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 1838baddfdSBarry Smith { 19dfbe8321SBarry Smith PetscErrorCode ierr; 2025023028SHong Zhang PetscBool scalable=PETSC_FALSE; 215c66b693SKris Buschelman 225c66b693SKris Buschelman PetscFunctionBegin; 2326be0446SHong Zhang if (scall == MAT_INITIAL_MATRIX){ 24598bc09dSHong Zhang ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr); 2525023028SHong Zhang ierr = PetscOptionsGetBool(PETSC_NULL,"-matmatmult_scalable",&scalable,PETSC_NULL);CHKERRQ(ierr); 2625023028SHong Zhang if (scalable){ 2725023028SHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(A,B,fill,C);CHKERRQ(ierr); 2825023028SHong Zhang } else { 2926be0446SHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr); 3025023028SHong Zhang } 31598bc09dSHong Zhang ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr); 3226be0446SHong Zhang } 33598bc09dSHong Zhang ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr); 3401e47db4SHong Zhang ierr = (*(*C)->ops->matmultnumeric)(A,B,*C);CHKERRQ(ierr); 35598bc09dSHong Zhang ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr); 365c66b693SKris Buschelman PetscFunctionReturn(0); 375c66b693SKris Buschelman } 38ec01deb9SMatthew Knepley EXTERN_C_END 391c24bd37SHong Zhang 40e9e4536cSHong Zhang /* 41e9e4536cSHong Zhang MatGetSymbolicMatMatMult_SeqAIJ_SeqAIJ - Get symbolic structure of C=A*B 42e9e4536cSHong Zhang Input Parameter: 43e9e4536cSHong Zhang . am, Ai, Aj - number of rows and structure of A 44e9e4536cSHong Zhang . bm, bn, Bi, Bj - number of rows, columns, and structure of B 45e9e4536cSHong Zhang . fill - filll ratio See MatMatMult() 46e9e4536cSHong Zhang 47e9e4536cSHong Zhang Output Parameter: 48e9e4536cSHong Zhang . Ci, Cj - structure of C = A*B 49e9e4536cSHong Zhang . nspacedouble - number of extra mallocs 50e9e4536cSHong Zhang */ 5126be0446SHong Zhang #undef __FUNCT__ 52b561aa9dSHong Zhang #define __FUNCT__ "MatGetSymbolicMatMatMult_SeqAIJ_SeqAIJ" 53b561aa9dSHong Zhang PetscErrorCode MatGetSymbolicMatMatMult_SeqAIJ_SeqAIJ(PetscInt am,PetscInt *Ai,PetscInt *Aj,PetscInt bm,PetscInt bn,PetscInt *Bi,PetscInt *Bj,PetscReal fill,PetscInt *Ci[],PetscInt *Cj[],PetscInt *nspacedouble) 5458c24d83SHong Zhang { 55dfbe8321SBarry Smith PetscErrorCode ierr; 56a1a86e44SBarry Smith PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 57b561aa9dSHong Zhang PetscInt *ai=Ai,*aj=Aj,*bi=Bi,*bj=Bj,*bjj,*ci,*cj; 58b561aa9dSHong Zhang PetscInt i,j,anzi,brow,bnzj,cnzi,nlnk,*lnk,ndouble=0; 59be0fcf8dSHong Zhang PetscBT lnkbt; 6058c24d83SHong Zhang 6158c24d83SHong Zhang PetscFunctionBegin; 6258c24d83SHong Zhang /* Allocate ci array, arrays for fill computation and */ 6358c24d83SHong Zhang /* free space for accumulating nonzero column info */ 6438baddfdSBarry Smith ierr = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr); 6558c24d83SHong Zhang ci[0] = 0; 6658c24d83SHong Zhang 67be0fcf8dSHong Zhang /* create and initialize a linked list */ 68be0fcf8dSHong Zhang nlnk = bn+1; 69be0fcf8dSHong Zhang ierr = PetscLLCreate(bn,bn,nlnk,lnk,lnkbt);CHKERRQ(ierr); 7058c24d83SHong Zhang 71c5db241fSHong Zhang /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */ 72a1a86e44SBarry Smith ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr); 7358c24d83SHong Zhang current_space = free_space; 7458c24d83SHong Zhang 7558c24d83SHong Zhang /* Determine symbolic info for each row of the product: */ 7658c24d83SHong Zhang for (i=0; i<am; i++) { 7758c24d83SHong Zhang anzi = ai[i+1] - ai[i]; 7858c24d83SHong Zhang cnzi = 0; 79b561aa9dSHong Zhang aj = Aj + ai[i]; 8077584fe6SHong Zhang for (j=0; j<anzi; j++){ 81b561aa9dSHong Zhang brow = aj[j]; 8258c24d83SHong Zhang bnzj = bi[brow+1] - bi[brow]; 8358c24d83SHong Zhang bjj = bj + bi[brow]; 841c239cc6SHong Zhang /* add non-zero cols of B into the sorted linked list lnk */ 85dadf0e6bSHong Zhang ierr = PetscLLAddSorted(bnzj,bjj,bn,nlnk,lnk,lnkbt);CHKERRQ(ierr); 861c239cc6SHong Zhang cnzi += nlnk; 8758c24d83SHong Zhang } 8858c24d83SHong Zhang 8958c24d83SHong Zhang /* If free space is not available, make more free space */ 9058c24d83SHong Zhang /* Double the amount of total space in the list */ 9158c24d83SHong Zhang if (current_space->local_remaining<cnzi) { 924238b7adSHong Zhang ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 93b561aa9dSHong Zhang ndouble++; 9458c24d83SHong Zhang } 9558c24d83SHong Zhang 96c5db241fSHong Zhang /* Copy data into free space, then initialize lnk */ 97be0fcf8dSHong Zhang ierr = PetscLLClean(bn,bn,cnzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 98c5db241fSHong Zhang current_space->array += cnzi; 9958c24d83SHong Zhang current_space->local_used += cnzi; 10058c24d83SHong Zhang current_space->local_remaining -= cnzi; 10158c24d83SHong Zhang ci[i+1] = ci[i] + cnzi; 10258c24d83SHong Zhang } 10358c24d83SHong Zhang 10458c24d83SHong Zhang /* Column indices are in the list of free space */ 10558c24d83SHong Zhang /* Allocate space for cj, initialize cj, and */ 10658c24d83SHong Zhang /* destroy list of free space and other temporary array(s) */ 10738baddfdSBarry Smith ierr = PetscMalloc((ci[am]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr); 108a1a86e44SBarry Smith ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr); 109be0fcf8dSHong Zhang ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 11058c24d83SHong Zhang 111b561aa9dSHong Zhang *Ci = ci; 112b561aa9dSHong Zhang *Cj = cj; 113b561aa9dSHong Zhang *nspacedouble = ndouble; 114b561aa9dSHong Zhang PetscFunctionReturn(0); 115b561aa9dSHong Zhang } 116b561aa9dSHong Zhang 11731e66dd0SHong Zhang /* 11825023028SHong Zhang MatGetSymbolicMatMatMult_SeqAIJ_SeqAIJ_Scalable - same as MatGetSymbolicMatMatMult_SeqAIJ_SeqAIJ() 11931e66dd0SHong Zhang but replaces O(bn) array 'lnkbt' with a scalable array 'abj' of size Crmax. 12031e66dd0SHong Zhang */ 121b561aa9dSHong Zhang #undef __FUNCT__ 12225023028SHong Zhang #define __FUNCT__ "MatGetSymbolicMatMatMult_SeqAIJ_SeqAIJ_Scalable" 12325023028SHong Zhang PetscErrorCode MatGetSymbolicMatMatMult_SeqAIJ_SeqAIJ_Scalable(PetscInt am,PetscInt *Ai,PetscInt *Aj,PetscInt bm,PetscInt bn,PetscInt *Bi,PetscInt *Bj,PetscReal fill,PetscInt *Ci[],PetscInt *Cj[],PetscInt *nspacedouble) 124e9e4536cSHong Zhang { 125e9e4536cSHong Zhang PetscErrorCode ierr; 12631e66dd0SHong Zhang PetscInt *aj=Aj,*bi=Bi,*bj,*ci,*cj,rmax=0,*abj,*cj_tmp,nextabj; 127e9e4536cSHong Zhang PetscInt i,j,anzi,brow,bnzj,cnzi,k; 128e9e4536cSHong Zhang PetscBT bt; 129164eb823SHong Zhang #if defined(DEBUG_MATMATMULT) 130164eb823SHong Zhang PetscLogDouble t0,tf,time_rmax=0.0,time_sym=0.0,time_sort=0.0,t0_sort,tf_sort; 131164eb823SHong Zhang #endif 132e9e4536cSHong Zhang 133e9e4536cSHong Zhang PetscFunctionBegin; 13431e66dd0SHong Zhang /* Allocate ci array and bt */ 135e9e4536cSHong Zhang ierr = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr); 136e9e4536cSHong Zhang ci[0] = 0; 137e9e4536cSHong Zhang 138164eb823SHong Zhang /* Get ci and count rmax of C - Crmax <= max(Armax*Brmax, bn) */ 139164eb823SHong Zhang #if defined(DEBUG_MATMATMULT) 140164eb823SHong Zhang ierr = PetscGetTime(&t0);CHKERRQ(ierr); 141164eb823SHong Zhang #endif 142e9e4536cSHong Zhang ierr = PetscBTCreate(bn,bt);CHKERRQ(ierr); 143e9e4536cSHong Zhang ierr = PetscBTMemzero(bn,bt);CHKERRQ(ierr); 144e9e4536cSHong Zhang for (i=0; i<am; i++) { 14531e66dd0SHong Zhang anzi = Ai[i+1] - Ai[i]; 146e9e4536cSHong Zhang cnzi = 0; 14731e66dd0SHong Zhang aj = Aj + Ai[i]; 148e9e4536cSHong Zhang for (j=0; j<anzi; j++){ 149e9e4536cSHong Zhang brow = aj[j]; 150e9e4536cSHong Zhang bnzj = bi[brow+1] - bi[brow]; 15131e66dd0SHong Zhang bj = Bj + bi[brow]; 152e9e4536cSHong Zhang for (k=0; k<bnzj; k++){ 15331e66dd0SHong Zhang if (!PetscBTLookupSet(bt,bj[k])) cnzi++; /* new entry */ 154e9e4536cSHong Zhang } 155e9e4536cSHong Zhang } 15631e66dd0SHong Zhang ierr = PetscBTMemzero(bn,bt);CHKERRQ(ierr); 157e9e4536cSHong Zhang ci[i+1] = ci[i] + cnzi; 158e9e4536cSHong Zhang if (rmax < cnzi) rmax = cnzi; 159e9e4536cSHong Zhang } 160164eb823SHong Zhang #if defined(DEBUG_MATMATMULT) 161164eb823SHong Zhang ierr = PetscGetTime(&tf);CHKERRQ(ierr); 162164eb823SHong Zhang time_rmax += tf - t0; 163164eb823SHong Zhang ierr = PetscGetTime(&t0);CHKERRQ(ierr); 164164eb823SHong Zhang #endif 165e9e4536cSHong Zhang /* Allocate space for cj */ 166e9e4536cSHong Zhang ierr = PetscMalloc((ci[am]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr); 167e9e4536cSHong Zhang 16831e66dd0SHong Zhang /* Allocate a temp array for storing column indices of A*B */ 169e9e4536cSHong Zhang ierr = PetscMalloc((rmax+1)*sizeof(PetscInt),&abj);CHKERRQ(ierr); 170e9e4536cSHong Zhang 171e9e4536cSHong Zhang /* Determine cj */ 172e9e4536cSHong Zhang for (i=0; i<am; i++) { 17331e66dd0SHong Zhang anzi = Ai[i+1] - Ai[i]; 174e9e4536cSHong Zhang cnzi = 0; 175e9e4536cSHong Zhang nextabj = 0; 17631e66dd0SHong Zhang aj = Aj + Ai[i]; 177e9e4536cSHong Zhang for (j=0; j<anzi; j++){ 178e9e4536cSHong Zhang brow = aj[j]; 179e9e4536cSHong Zhang bnzj = bi[brow+1] - bi[brow]; 18031e66dd0SHong Zhang bj = Bj + bi[brow]; 181e9e4536cSHong Zhang for (k=0; k<bnzj; k++){ 18231e66dd0SHong Zhang if (!PetscBTLookupSet(bt,bj[k])){ /* new entry */ 18331e66dd0SHong Zhang abj[nextabj] = bj[k]; nextabj++; 184e9e4536cSHong Zhang } 185e9e4536cSHong Zhang } 186e9e4536cSHong Zhang } 187e9e4536cSHong Zhang 18831e66dd0SHong Zhang /* Sort abj, then copy it to cj */ 189e9e4536cSHong Zhang cnzi = ci[i+1] - ci[i]; 190164eb823SHong Zhang #if defined(DEBUG_MATMATMULT) 191164eb823SHong Zhang ierr = PetscGetTime(&t0_sort);CHKERRQ(ierr); 192164eb823SHong Zhang #endif 193e9e4536cSHong Zhang ierr = PetscSortInt(cnzi,abj);CHKERRQ(ierr); 194164eb823SHong Zhang #if defined(DEBUG_MATMATMULT) 195164eb823SHong Zhang ierr = PetscGetTime(&tf_sort);CHKERRQ(ierr); 196164eb823SHong Zhang time_sort += tf_sort - t0_sort; 197164eb823SHong Zhang #endif 198164eb823SHong Zhang 199e9e4536cSHong Zhang cj_tmp = cj + ci[i]; 200e9e4536cSHong Zhang for (k=0; k< cnzi; k++){ 201e9e4536cSHong Zhang cj_tmp[k] = abj[k]; 20231e66dd0SHong Zhang ierr = PetscBTClear(bt,abj[k]);CHKERRQ(ierr); 203e9e4536cSHong Zhang } 204e9e4536cSHong Zhang } 205164eb823SHong Zhang #if defined(DEBUG_MATMATMULT) 206164eb823SHong Zhang ierr = PetscGetTime(&tf);CHKERRQ(ierr); 207164eb823SHong Zhang time_sym += tf - t0; 208164eb823SHong Zhang #endif 209e9e4536cSHong Zhang 210e9e4536cSHong Zhang ierr = PetscBTDestroy(bt);CHKERRQ(ierr); 211e9e4536cSHong Zhang ierr = PetscFree(abj);CHKERRQ(ierr); 212e9e4536cSHong Zhang 213e9e4536cSHong Zhang *Ci = ci; 214e9e4536cSHong Zhang *Cj = cj; 215e9e4536cSHong Zhang *nspacedouble = 0; 216164eb823SHong Zhang #if defined(DEBUG_MATMATMULT) 217164eb823SHong Zhang printf("Time of GetSymbolic_Scalable: rmax %g + sym %g (sort %g) = %g; Crmax %d, pN %d\n",time_rmax,time_sym,time_sort,time_rmax+time_sym,rmax,bn); 218164eb823SHong Zhang #endif 219e9e4536cSHong Zhang PetscFunctionReturn(0); 220e9e4536cSHong Zhang } 221e9e4536cSHong Zhang 222e9e4536cSHong Zhang #undef __FUNCT__ 223b561aa9dSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ" 224b561aa9dSHong Zhang PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C) 225b561aa9dSHong Zhang { 226b561aa9dSHong Zhang PetscErrorCode ierr; 227b561aa9dSHong Zhang Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c; 228b561aa9dSHong Zhang PetscInt *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*ci,*cj; 229b561aa9dSHong Zhang PetscInt am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N,nspacedouble; 230b561aa9dSHong Zhang MatScalar *ca; 231b561aa9dSHong Zhang PetscReal afill; 232b561aa9dSHong Zhang 233b561aa9dSHong Zhang PetscFunctionBegin; 234b561aa9dSHong Zhang /* Get ci and cj */ 235b561aa9dSHong Zhang ierr = MatGetSymbolicMatMatMult_SeqAIJ_SeqAIJ(am,ai,aj,bm,bn,bi,bj,fill,&ci,&cj,&nspacedouble);CHKERRQ(ierr); 236e9e4536cSHong Zhang #if defined(DEBUG_MATMATMULT) 237e9e4536cSHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF,"MatGetSymbolicMatMatMult_SeqAIJ_SeqAIJ() is done \n");CHKERRQ(ierr); 238e9e4536cSHong Zhang #endif 239b561aa9dSHong Zhang 24058c24d83SHong Zhang /* Allocate space for ca */ 24158c24d83SHong Zhang ierr = PetscMalloc((ci[am]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr); 24258c24d83SHong Zhang ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr); 24358c24d83SHong Zhang 24426be0446SHong Zhang /* put together the new symbolic matrix */ 2457adad957SLisandro Dalcin ierr = MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,am,bn,ci,cj,ca,C);CHKERRQ(ierr); 24658c24d83SHong Zhang 24758c24d83SHong Zhang /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 24858c24d83SHong Zhang /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */ 24958c24d83SHong Zhang c = (Mat_SeqAIJ *)((*C)->data); 250e6b907acSBarry Smith c->free_a = PETSC_TRUE; 251e6b907acSBarry Smith c->free_ij = PETSC_TRUE; 25258c24d83SHong Zhang c->nonew = 0; 2538cdbd757SHong Zhang (*C)->ops->matmult = MatMatMult_SeqAIJ_SeqAIJ; 25458c24d83SHong Zhang 25536ec6d2dSHong Zhang ierr = PetscMalloc(bn*sizeof(PetscScalar),&c->matmult_abdense);CHKERRQ(ierr); 25625023028SHong Zhang ierr = PetscMemzero(c->matmult_abdense,bn*sizeof(PetscScalar));CHKERRQ(ierr); 257c58ca2e3SHong Zhang (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ; /* fast, takes additional dense_axpy*bn*sizeof(PetscScalar) space */ 2580b7e3e3dSHong Zhang 2597212da7cSHong Zhang /* set MatInfo */ 2607212da7cSHong Zhang afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5; 261f2b054eeSHong Zhang if (afill < 1.0) afill = 1.0; 2627212da7cSHong Zhang c->maxnz = ci[am]; 2637212da7cSHong Zhang c->nz = ci[am]; 2647212da7cSHong Zhang (*C)->info.mallocs = nspacedouble; 2657212da7cSHong Zhang (*C)->info.fill_ratio_given = fill; 2667212da7cSHong Zhang (*C)->info.fill_ratio_needed = afill; 2677212da7cSHong Zhang 2687212da7cSHong Zhang #if defined(PETSC_USE_INFO) 2697212da7cSHong Zhang if (ci[am]) { 270f2b054eeSHong Zhang ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr); 271f2b054eeSHong Zhang ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr); 272f2b054eeSHong Zhang } else { 273f2b054eeSHong Zhang ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr); 274be0fcf8dSHong Zhang } 275f2b054eeSHong Zhang #endif 27658c24d83SHong Zhang PetscFunctionReturn(0); 27758c24d83SHong Zhang } 278d50806bdSBarry Smith 27926be0446SHong Zhang #undef __FUNCT__ 28026be0446SHong Zhang #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqAIJ" 281dfbe8321SBarry Smith PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C) 282d50806bdSBarry Smith { 283dfbe8321SBarry Smith PetscErrorCode ierr; 284d13dce4bSSatish Balay PetscLogDouble flops=0.0; 285d50806bdSBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 286d50806bdSBarry Smith Mat_SeqAIJ *b = (Mat_SeqAIJ *)B->data; 287d50806bdSBarry Smith Mat_SeqAIJ *c = (Mat_SeqAIJ *)C->data; 28838baddfdSBarry Smith PetscInt *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j; 289c58ca2e3SHong Zhang PetscInt am=A->rmap->n,cm=C->rmap->n; 290519eb980SHong Zhang PetscInt i,j,k,anzi,bnzi,cnzi,brow; 29136ec6d2dSHong Zhang PetscScalar *aa=a->a,*ba=b->a,*baj,*ca=c->a,valtmp; 2920b7e3e3dSHong Zhang PetscScalar *ab_dense=c->matmult_abdense; 293d50806bdSBarry Smith 294d50806bdSBarry Smith PetscFunctionBegin; 295c124e916SHong Zhang /* clean old values in C */ 296c124e916SHong Zhang ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr); 297d50806bdSBarry Smith /* Traverse A row-wise. */ 298d50806bdSBarry Smith /* Build the ith row in C by summing over nonzero columns in A, */ 299d50806bdSBarry Smith /* the rows of B corresponding to nonzeros of A. */ 300d50806bdSBarry Smith for (i=0; i<am; i++) { 301d50806bdSBarry Smith anzi = ai[i+1] - ai[i]; 302d50806bdSBarry Smith for (j=0; j<anzi; j++) { 303519eb980SHong Zhang brow = aj[j]; 304d50806bdSBarry Smith bnzi = bi[brow+1] - bi[brow]; 305d50806bdSBarry Smith bjj = bj + bi[brow]; 306d50806bdSBarry Smith baj = ba + bi[brow]; 307519eb980SHong Zhang /* perform dense axpy */ 30836ec6d2dSHong Zhang valtmp = aa[j]; 309519eb980SHong Zhang for (k=0; k<bnzi; k++) { 31036ec6d2dSHong Zhang ab_dense[bjj[k]] += valtmp*baj[k]; 311519eb980SHong Zhang } 312519eb980SHong Zhang flops += 2*bnzi; 313519eb980SHong Zhang } 314c58ca2e3SHong Zhang aj += anzi; aa += anzi; 315c58ca2e3SHong Zhang 316c58ca2e3SHong Zhang cnzi = ci[i+1] - ci[i]; 317519eb980SHong Zhang for (k=0; k<cnzi; k++) { 318519eb980SHong Zhang ca[k] += ab_dense[cj[k]]; 319519eb980SHong Zhang ab_dense[cj[k]] = 0.0; /* zero ab_dense */ 320519eb980SHong Zhang } 321519eb980SHong Zhang flops += cnzi; 322519eb980SHong Zhang cj += cnzi; ca += cnzi; 323519eb980SHong Zhang } 324c58ca2e3SHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 325c58ca2e3SHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 326c58ca2e3SHong Zhang ierr = PetscLogFlops(flops);CHKERRQ(ierr); 327c58ca2e3SHong Zhang PetscFunctionReturn(0); 328c58ca2e3SHong Zhang } 329c58ca2e3SHong Zhang 330c58ca2e3SHong Zhang #undef __FUNCT__ 33125023028SHong Zhang #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable" 33225023028SHong Zhang PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(Mat A,Mat B,Mat C) 333c58ca2e3SHong Zhang { 334c58ca2e3SHong Zhang PetscErrorCode ierr; 335c58ca2e3SHong Zhang PetscLogDouble flops=0.0; 336c58ca2e3SHong Zhang Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 337c58ca2e3SHong Zhang Mat_SeqAIJ *b = (Mat_SeqAIJ *)B->data; 338c58ca2e3SHong Zhang Mat_SeqAIJ *c = (Mat_SeqAIJ *)C->data; 339c58ca2e3SHong Zhang PetscInt *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j; 340c58ca2e3SHong Zhang PetscInt am=A->rmap->N,cm=C->rmap->N; 341c58ca2e3SHong Zhang PetscInt i,j,k,anzi,bnzi,cnzi,brow; 34236ec6d2dSHong Zhang PetscScalar *aa=a->a,*ba=b->a,*baj,*ca=c->a,valtmp; 343c58ca2e3SHong Zhang PetscInt nextb; 344c58ca2e3SHong Zhang 345c58ca2e3SHong Zhang PetscFunctionBegin; 346e9e4536cSHong Zhang #if defined(DEBUG_MATMATMULT) 34725023028SHong Zhang //ierr = PetscPrintf(PETSC_COMM_SELF,"MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable...\n");CHKERRQ(ierr); 348e9e4536cSHong Zhang #endif 349c58ca2e3SHong Zhang /* clean old values in C */ 350c58ca2e3SHong Zhang ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr); 351c58ca2e3SHong Zhang /* Traverse A row-wise. */ 352c58ca2e3SHong Zhang /* Build the ith row in C by summing over nonzero columns in A, */ 353c58ca2e3SHong Zhang /* the rows of B corresponding to nonzeros of A. */ 354519eb980SHong Zhang for (i=0;i<am;i++) { 355519eb980SHong Zhang anzi = ai[i+1] - ai[i]; 356519eb980SHong Zhang cnzi = ci[i+1] - ci[i]; 357519eb980SHong Zhang for (j=0;j<anzi;j++) { 358519eb980SHong Zhang brow = aj[j]; 359519eb980SHong Zhang bnzi = bi[brow+1] - bi[brow]; 360519eb980SHong Zhang bjj = bj + bi[brow]; 361519eb980SHong Zhang baj = ba + bi[brow]; 362519eb980SHong Zhang /* perform sparse axpy */ 36336ec6d2dSHong Zhang valtmp = aa[j]; 364c124e916SHong Zhang nextb = 0; 365c124e916SHong Zhang for (k=0; nextb<bnzi; k++) { 366c124e916SHong Zhang if (cj[k] == bjj[nextb]){ /* ccol == bcol */ 36736ec6d2dSHong Zhang ca[k] += valtmp*baj[nextb++]; 368c124e916SHong Zhang } 369d50806bdSBarry Smith } 370d50806bdSBarry Smith flops += 2*bnzi; 371d50806bdSBarry Smith } 372519eb980SHong Zhang aj += anzi; aa += anzi; 373519eb980SHong Zhang cj += cnzi; ca += cnzi; 374d50806bdSBarry Smith } 375c58ca2e3SHong Zhang 376716bacf3SKris Buschelman ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 377716bacf3SKris Buschelman ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 378d50806bdSBarry Smith ierr = PetscLogFlops(flops);CHKERRQ(ierr); 379d50806bdSBarry Smith PetscFunctionReturn(0); 380d50806bdSBarry Smith } 381bc011b1eSHong Zhang 382e9e4536cSHong Zhang #undef __FUNCT__ 38325023028SHong Zhang #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable" 38425023028SHong Zhang PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(Mat A,Mat B,PetscReal fill,Mat *C) 385e9e4536cSHong Zhang { 386e9e4536cSHong Zhang PetscErrorCode ierr; 387e9e4536cSHong Zhang Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c; 388e9e4536cSHong Zhang PetscInt *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*ci,*cj; 389e9e4536cSHong Zhang PetscInt am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N,nspacedouble; 390e9e4536cSHong Zhang MatScalar *ca; 391e9e4536cSHong Zhang PetscReal afill; 392e9e4536cSHong Zhang 393e9e4536cSHong Zhang PetscFunctionBegin; 394e9e4536cSHong Zhang #if defined(DEBUG_MATMATMULT) 395164eb823SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF,"MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable Armax %d, Brmax %d\n",a->rmax,b->rmax);CHKERRQ(ierr); 396e9e4536cSHong Zhang #endif 397e9e4536cSHong Zhang /* Get ci and cj */ 39825023028SHong Zhang ierr = MatGetSymbolicMatMatMult_SeqAIJ_SeqAIJ_Scalable(am,ai,aj,bm,bn,bi,bj,fill,&ci,&cj,&nspacedouble);CHKERRQ(ierr); 399e9e4536cSHong Zhang #if defined(DEBUG_MATMATMULT) 40025023028SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF,"MatGetSymbolicMatMatMult_SeqAIJ_SeqAIJ_Scalable() is done \n");CHKERRQ(ierr); 401e9e4536cSHong Zhang #endif 402e9e4536cSHong Zhang 403e9e4536cSHong Zhang /* Allocate space for ca */ 404e9e4536cSHong Zhang ierr = PetscMalloc((ci[am]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr); 405e9e4536cSHong Zhang ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr); 406e9e4536cSHong Zhang 407e9e4536cSHong Zhang /* put together the new symbolic matrix */ 408e9e4536cSHong Zhang ierr = MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,am,bn,ci,cj,ca,C);CHKERRQ(ierr); 409e9e4536cSHong Zhang 410e9e4536cSHong Zhang /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 411e9e4536cSHong Zhang /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */ 412e9e4536cSHong Zhang c = (Mat_SeqAIJ *)((*C)->data); 413e9e4536cSHong Zhang c->free_a = PETSC_TRUE; 414e9e4536cSHong Zhang c->free_ij = PETSC_TRUE; 415e9e4536cSHong Zhang c->nonew = 0; 41625023028SHong Zhang (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable; /* slower, less memory */ 417e9e4536cSHong Zhang 418e9e4536cSHong Zhang /* set MatInfo */ 419e9e4536cSHong Zhang afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5; 420e9e4536cSHong Zhang if (afill < 1.0) afill = 1.0; 421e9e4536cSHong Zhang c->maxnz = ci[am]; 422e9e4536cSHong Zhang c->nz = ci[am]; 423e9e4536cSHong Zhang (*C)->info.mallocs = nspacedouble; 424e9e4536cSHong Zhang (*C)->info.fill_ratio_given = fill; 425e9e4536cSHong Zhang (*C)->info.fill_ratio_needed = afill; 426e9e4536cSHong Zhang 427e9e4536cSHong Zhang #if defined(PETSC_USE_INFO) 428e9e4536cSHong Zhang if (ci[am]) { 429e9e4536cSHong Zhang ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr); 430e9e4536cSHong Zhang ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr); 431e9e4536cSHong Zhang } else { 432e9e4536cSHong Zhang ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr); 433e9e4536cSHong Zhang } 434e9e4536cSHong Zhang #endif 435e9e4536cSHong Zhang PetscFunctionReturn(0); 436e9e4536cSHong Zhang } 437e9e4536cSHong Zhang 438e9e4536cSHong Zhang 439d2853540SHong Zhang /* This routine is not used. Should be removed! */ 440bc011b1eSHong Zhang #undef __FUNCT__ 4416fc122caSHong Zhang #define __FUNCT__ "MatMatTransposeMult_SeqAIJ_SeqAIJ" 4426fc122caSHong Zhang PetscErrorCode MatMatTransposeMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 4435df89d91SHong Zhang { 444bc011b1eSHong Zhang PetscErrorCode ierr; 445bc011b1eSHong Zhang 446bc011b1eSHong Zhang PetscFunctionBegin; 447bc011b1eSHong Zhang if (scall == MAT_INITIAL_MATRIX){ 4486fc122caSHong Zhang ierr = MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr); 449bc011b1eSHong Zhang } 4506fc122caSHong Zhang ierr = MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr); 451bc011b1eSHong Zhang PetscFunctionReturn(0); 452bc011b1eSHong Zhang } 453bc011b1eSHong Zhang 454bc011b1eSHong Zhang #undef __FUNCT__ 455e286af10SHong Zhang #define __FUNCT__ "PetscContainerDestroy_Mat_MatMatTransMult" 456e286af10SHong Zhang PetscErrorCode PetscContainerDestroy_Mat_MatMatTransMult(void *ptr) 4572128a86cSHong Zhang { 4582128a86cSHong Zhang PetscErrorCode ierr; 459e286af10SHong Zhang Mat_MatMatTransMult *multtrans=(Mat_MatMatTransMult*)ptr; 4602128a86cSHong Zhang 4612128a86cSHong Zhang PetscFunctionBegin; 462b9af6bddSHong Zhang ierr = MatTransposeColoringDestroy(&multtrans->matcoloring);CHKERRQ(ierr); 4632128a86cSHong Zhang ierr = MatDestroy(&multtrans->Bt_den);CHKERRQ(ierr); 4642128a86cSHong Zhang ierr = MatDestroy(&multtrans->ABt_den);CHKERRQ(ierr); 4652128a86cSHong Zhang ierr = PetscFree(multtrans);CHKERRQ(ierr); 4662128a86cSHong Zhang PetscFunctionReturn(0); 4672128a86cSHong Zhang } 4682128a86cSHong Zhang 4692128a86cSHong Zhang #undef __FUNCT__ 4702128a86cSHong Zhang #define __FUNCT__ "MatDestroy_SeqAIJ_MatMatMultTrans" 4712128a86cSHong Zhang PetscErrorCode MatDestroy_SeqAIJ_MatMatMultTrans(Mat A) 4722128a86cSHong Zhang { 4732128a86cSHong Zhang PetscErrorCode ierr; 4742128a86cSHong Zhang PetscContainer container; 475e286af10SHong Zhang Mat_MatMatTransMult *multtrans=PETSC_NULL; 4762128a86cSHong Zhang 4772128a86cSHong Zhang PetscFunctionBegin; 478e286af10SHong Zhang ierr = PetscObjectQuery((PetscObject)A,"Mat_MatMatTransMult",(PetscObject *)&container);CHKERRQ(ierr); 4792128a86cSHong Zhang if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit"); 4802128a86cSHong Zhang ierr = PetscContainerGetPointer(container,(void **)&multtrans);CHKERRQ(ierr); 4812128a86cSHong Zhang A->ops->destroy = multtrans->destroy; 4822128a86cSHong Zhang if (A->ops->destroy) { 4832128a86cSHong Zhang ierr = (*A->ops->destroy)(A);CHKERRQ(ierr); 4842128a86cSHong Zhang } 485e286af10SHong Zhang ierr = PetscObjectCompose((PetscObject)A,"Mat_MatMatTransMult",0);CHKERRQ(ierr); 4862128a86cSHong Zhang PetscFunctionReturn(0); 4872128a86cSHong Zhang } 4882128a86cSHong Zhang 4892128a86cSHong Zhang #undef __FUNCT__ 4906fc122caSHong Zhang #define __FUNCT__ "MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ" 4916fc122caSHong Zhang PetscErrorCode MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C) 492bc011b1eSHong Zhang { 4935df89d91SHong Zhang PetscErrorCode ierr; 49481d82fe4SHong Zhang Mat Bt; 49581d82fe4SHong Zhang PetscInt *bti,*btj; 496e286af10SHong Zhang Mat_MatMatTransMult *multtrans; 4972128a86cSHong Zhang PetscContainer container; 498d2853540SHong Zhang PetscLogDouble t0,tf,etime2=0.0; 499d2853540SHong Zhang 50081d82fe4SHong Zhang PetscFunctionBegin; 501d2853540SHong Zhang ierr = PetscGetTime(&t0);CHKERRQ(ierr); 50281d82fe4SHong Zhang /* create symbolic Bt */ 50381d82fe4SHong Zhang ierr = MatGetSymbolicTranspose_SeqAIJ(B,&bti,&btj);CHKERRQ(ierr); 50481d82fe4SHong Zhang ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,B->cmap->n,B->rmap->n,bti,btj,PETSC_NULL,&Bt);CHKERRQ(ierr); 50581d82fe4SHong Zhang 50681d82fe4SHong Zhang /* get symbolic C=A*Bt */ 50781d82fe4SHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,Bt,fill,C);CHKERRQ(ierr); 50881d82fe4SHong Zhang 5092128a86cSHong Zhang /* create a supporting struct for reuse intermidiate dense matrices with matcoloring */ 510e286af10SHong Zhang ierr = PetscNew(Mat_MatMatTransMult,&multtrans);CHKERRQ(ierr); 5112128a86cSHong Zhang 5122128a86cSHong Zhang /* attach the supporting struct to C */ 5132128a86cSHong Zhang ierr = PetscContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr); 5142128a86cSHong Zhang ierr = PetscContainerSetPointer(container,multtrans);CHKERRQ(ierr); 515e286af10SHong Zhang ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_Mat_MatMatTransMult);CHKERRQ(ierr); 516e286af10SHong Zhang ierr = PetscObjectCompose((PetscObject)(*C),"Mat_MatMatTransMult",(PetscObject)container);CHKERRQ(ierr); 5172128a86cSHong Zhang ierr = PetscContainerDestroy(&container);CHKERRQ(ierr); 5182128a86cSHong Zhang 5192128a86cSHong Zhang multtrans->usecoloring = PETSC_FALSE; 5202128a86cSHong Zhang multtrans->destroy = (*C)->ops->destroy; 5212128a86cSHong Zhang (*C)->ops->destroy = MatDestroy_SeqAIJ_MatMatMultTrans; 5222128a86cSHong Zhang 523d2853540SHong Zhang ierr = PetscGetTime(&tf);CHKERRQ(ierr); 524d2853540SHong Zhang etime2 += tf - t0; 525d2853540SHong Zhang 526f400d4cbSHong Zhang ierr = PetscOptionsGetBool(PETSC_NULL,"-matmattransmult_color",&multtrans->usecoloring,PETSC_NULL);CHKERRQ(ierr); 5272128a86cSHong Zhang if (multtrans->usecoloring){ 528b9af6bddSHong Zhang /* Create MatTransposeColoring from symbolic C=A*B^T */ 529b9af6bddSHong Zhang MatTransposeColoring matcoloring; 5302128a86cSHong Zhang ISColoring iscoloring; 5312128a86cSHong Zhang Mat Bt_dense,C_dense; 532d2853540SHong Zhang PetscLogDouble etime0=0.0,etime01=0.0,etime1=0.0; 533e8354b3eSHong Zhang 534e8354b3eSHong Zhang ierr = PetscGetTime(&t0);CHKERRQ(ierr); 535502de53fSHong Zhang ierr = MatGetColoring(*C,MATCOLORINGLF,&iscoloring);CHKERRQ(ierr); 536d2853540SHong Zhang ierr = PetscGetTime(&tf);CHKERRQ(ierr); 537d2853540SHong Zhang etime0 += tf - t0; 538d2853540SHong Zhang 539d2853540SHong Zhang ierr = PetscGetTime(&t0);CHKERRQ(ierr); 540b9af6bddSHong Zhang ierr = MatTransposeColoringCreate(*C,iscoloring,&matcoloring);CHKERRQ(ierr); 5412128a86cSHong Zhang multtrans->matcoloring = matcoloring; 5422128a86cSHong Zhang ierr = ISColoringDestroy(&iscoloring);CHKERRQ(ierr); 543e8354b3eSHong Zhang ierr = PetscGetTime(&tf);CHKERRQ(ierr); 544d2853540SHong Zhang etime01 += tf - t0; 5452128a86cSHong Zhang 546e8354b3eSHong Zhang ierr = PetscGetTime(&t0);CHKERRQ(ierr); 5472128a86cSHong Zhang /* Create Bt_dense and C_dense = A*Bt_dense */ 5482128a86cSHong Zhang ierr = MatCreate(PETSC_COMM_SELF,&Bt_dense);CHKERRQ(ierr); 5492128a86cSHong Zhang ierr = MatSetSizes(Bt_dense,A->cmap->n,matcoloring->ncolors,A->cmap->n,matcoloring->ncolors);CHKERRQ(ierr); 5502128a86cSHong Zhang ierr = MatSetType(Bt_dense,MATSEQDENSE);CHKERRQ(ierr); 5512128a86cSHong Zhang ierr = MatSeqDenseSetPreallocation(Bt_dense,PETSC_NULL);CHKERRQ(ierr); 5522128a86cSHong Zhang Bt_dense->assembled = PETSC_TRUE; 5532128a86cSHong Zhang multtrans->Bt_den = Bt_dense; 5542128a86cSHong Zhang 5552128a86cSHong Zhang ierr = MatCreate(PETSC_COMM_SELF,&C_dense);CHKERRQ(ierr); 5562128a86cSHong Zhang ierr = MatSetSizes(C_dense,A->rmap->n,matcoloring->ncolors,A->rmap->n,matcoloring->ncolors);CHKERRQ(ierr); 5572128a86cSHong Zhang ierr = MatSetType(C_dense,MATSEQDENSE);CHKERRQ(ierr); 5582128a86cSHong Zhang ierr = MatSeqDenseSetPreallocation(C_dense,PETSC_NULL);CHKERRQ(ierr); 5592128a86cSHong Zhang Bt_dense->assembled = PETSC_TRUE; 5602128a86cSHong Zhang multtrans->ABt_den = C_dense; 561e8354b3eSHong Zhang ierr = PetscGetTime(&tf);CHKERRQ(ierr); 562e8354b3eSHong Zhang etime1 += tf - t0; 563f94ccd6cSHong Zhang 564f94ccd6cSHong Zhang #if defined(PETSC_USE_INFO) 5651ce71dffSSatish Balay { 566f94ccd6cSHong Zhang Mat_SeqAIJ *c=(Mat_SeqAIJ*)(*C)->data; 567f94ccd6cSHong 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)); 568f94ccd6cSHong Zhang ierr = PetscInfo5(*C,"Sym = GetColor %g + ColorCreate %g + MatDenseCreate %g + non-colorSym %g = %g\n",etime0,etime01,etime1,etime2,etime0+etime01+etime1+etime2); 5691ce71dffSSatish Balay } 570f94ccd6cSHong Zhang #endif 5712128a86cSHong Zhang } 572*e99be685SHong Zhang /* clean up */ 573*e99be685SHong Zhang ierr = MatDestroy(&Bt);CHKERRQ(ierr); 574*e99be685SHong Zhang ierr = MatRestoreSymbolicTranspose_SeqAIJ(B,&bti,&btj);CHKERRQ(ierr); 5752128a86cSHong Zhang 576f94ccd6cSHong Zhang 577f94ccd6cSHong Zhang 57881d82fe4SHong Zhang #if defined(INEFFICIENT_ALGORITHM) 57981d82fe4SHong Zhang /* The algorithm below computes am*bm sparse inner-product - inefficient! It will be deleted later. */ 5801fa1209cSHong Zhang PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 5811fa1209cSHong Zhang Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c; 5821fa1209cSHong Zhang PetscInt *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*ci,*cj,*acol,*bcol; 5831fa1209cSHong Zhang PetscInt am=A->rmap->N,bm=B->rmap->N; 5841fa1209cSHong Zhang PetscInt i,j,anzi,bnzj,cnzi,nlnk,*lnk,nspacedouble=0,ka,kb,index[1]; 5851fa1209cSHong Zhang MatScalar *ca; 5861fa1209cSHong Zhang PetscBT lnkbt; 5871fa1209cSHong Zhang PetscReal afill; 5885df89d91SHong Zhang 5891fa1209cSHong Zhang /* Allocate row pointer array ci */ 5901fa1209cSHong Zhang ierr = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr); 5911fa1209cSHong Zhang ci[0] = 0; 5921fa1209cSHong Zhang 5931fa1209cSHong Zhang /* Create and initialize a linked list for C columns */ 5941fa1209cSHong Zhang nlnk = bm+1; 5951fa1209cSHong Zhang ierr = PetscLLCreate(bm,bm,nlnk,lnk,lnkbt);CHKERRQ(ierr); 5961fa1209cSHong Zhang 5971fa1209cSHong Zhang /* Initial FreeSpace with size fill*(nnz(A)+nnz(B)) */ 5981fa1209cSHong Zhang ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr); 5991fa1209cSHong Zhang current_space = free_space; 6001fa1209cSHong Zhang 6011fa1209cSHong Zhang /* Determine symbolic info for each row of the product A*B^T: */ 6021fa1209cSHong Zhang for (i=0; i<am; i++) { 6031fa1209cSHong Zhang anzi = ai[i+1] - ai[i]; 6041fa1209cSHong Zhang cnzi = 0; 6051fa1209cSHong Zhang acol = aj + ai[i]; 6061fa1209cSHong Zhang for (j=0; j<bm; j++){ 6071fa1209cSHong Zhang bnzj = bi[j+1] - bi[j]; 6081fa1209cSHong Zhang bcol= bj + bi[j]; 60981d82fe4SHong Zhang /* sparse inner-product c(i,j)=A[i,:]*B[j,:]^T */ 6101fa1209cSHong Zhang ka = 0; kb = 0; 6111fa1209cSHong Zhang while (ka < anzi && kb < bnzj){ 61281d82fe4SHong Zhang while (acol[ka] < bcol[kb] && ka < anzi) ka++; 61381d82fe4SHong Zhang if (ka == anzi) break; 61481d82fe4SHong Zhang while (acol[ka] > bcol[kb] && kb < bnzj) kb++; 61581d82fe4SHong Zhang if (kb == bnzj) break; 6161fa1209cSHong Zhang if (acol[ka] == bcol[kb]){ /* add nonzero c(i,j) to lnk */ 6171fa1209cSHong Zhang index[0] = j; 6181fa1209cSHong Zhang ierr = PetscLLAdd(1,index,bm,nlnk,lnk,lnkbt);CHKERRQ(ierr); 6191fa1209cSHong Zhang cnzi++; 6201fa1209cSHong Zhang break; 6211fa1209cSHong Zhang } 6221fa1209cSHong Zhang } 6231fa1209cSHong Zhang } 6241fa1209cSHong Zhang 6251fa1209cSHong Zhang /* If free space is not available, make more free space */ 6261fa1209cSHong Zhang /* Double the amount of total space in the list */ 6271fa1209cSHong Zhang if (current_space->local_remaining<cnzi) { 6281fa1209cSHong Zhang ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 6291fa1209cSHong Zhang nspacedouble++; 6301fa1209cSHong Zhang } 6311fa1209cSHong Zhang 6321fa1209cSHong Zhang /* Copy data into free space, then initialize lnk */ 6331fa1209cSHong Zhang ierr = PetscLLClean(bm,bm,cnzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 6341fa1209cSHong Zhang current_space->array += cnzi; 6351fa1209cSHong Zhang current_space->local_used += cnzi; 6361fa1209cSHong Zhang current_space->local_remaining -= cnzi; 6371fa1209cSHong Zhang 6381fa1209cSHong Zhang ci[i+1] = ci[i] + cnzi; 6391fa1209cSHong Zhang } 6401fa1209cSHong Zhang 6411fa1209cSHong Zhang 6421fa1209cSHong Zhang /* Column indices are in the list of free space. 6431fa1209cSHong Zhang Allocate array cj, copy column indices to cj, and destroy list of free space */ 6441fa1209cSHong Zhang ierr = PetscMalloc((ci[am]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr); 6451fa1209cSHong Zhang ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr); 6461fa1209cSHong Zhang ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 6471fa1209cSHong Zhang 6481fa1209cSHong Zhang /* Allocate space for ca */ 6491fa1209cSHong Zhang ierr = PetscMalloc((ci[am]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr); 6501fa1209cSHong Zhang ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr); 6511fa1209cSHong Zhang 6521fa1209cSHong Zhang /* put together the new symbolic matrix */ 6531fa1209cSHong Zhang ierr = MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,am,bm,ci,cj,ca,C);CHKERRQ(ierr); 6541fa1209cSHong Zhang 6551fa1209cSHong Zhang /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 6561fa1209cSHong Zhang /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */ 6571fa1209cSHong Zhang c = (Mat_SeqAIJ *)((*C)->data); 6581fa1209cSHong Zhang c->free_a = PETSC_TRUE; 6591fa1209cSHong Zhang c->free_ij = PETSC_TRUE; 6601fa1209cSHong Zhang c->nonew = 0; 6611fa1209cSHong Zhang 6621fa1209cSHong Zhang /* set MatInfo */ 6631fa1209cSHong Zhang afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5; 6641fa1209cSHong Zhang if (afill < 1.0) afill = 1.0; 6651fa1209cSHong Zhang c->maxnz = ci[am]; 6661fa1209cSHong Zhang c->nz = ci[am]; 6671fa1209cSHong Zhang (*C)->info.mallocs = nspacedouble; 6681fa1209cSHong Zhang (*C)->info.fill_ratio_given = fill; 6691fa1209cSHong Zhang (*C)->info.fill_ratio_needed = afill; 6701fa1209cSHong Zhang 6711fa1209cSHong Zhang #if defined(PETSC_USE_INFO) 6721fa1209cSHong Zhang if (ci[am]) { 6731fa1209cSHong Zhang ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr); 6746fc122caSHong Zhang ierr = PetscInfo1((*C),"Use MatMatTransposeMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr); 6751fa1209cSHong Zhang } else { 6761fa1209cSHong Zhang ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr); 6771fa1209cSHong Zhang } 6781fa1209cSHong Zhang #endif 67981d82fe4SHong Zhang #endif 6805df89d91SHong Zhang PetscFunctionReturn(0); 6815df89d91SHong Zhang } 6825df89d91SHong Zhang 6836973769bSHong Zhang /* #define USE_ARRAY - for sparse dot product. Slower than !USE_ARRAY */ 6845df89d91SHong Zhang #undef __FUNCT__ 6856fc122caSHong Zhang #define __FUNCT__ "MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ" 6866fc122caSHong Zhang PetscErrorCode MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C) 6875df89d91SHong Zhang { 6885df89d91SHong Zhang PetscErrorCode ierr; 6895df89d91SHong Zhang Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data; 690e2cac8caSJed Brown PetscInt *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,anzi,bnzj,nexta,nextb,*acol,*bcol,brow; 69181d82fe4SHong Zhang PetscInt cm=C->rmap->n,*ci=c->i,*cj=c->j,i,j,cnzi,*ccol; 6925df89d91SHong Zhang PetscLogDouble flops=0.0; 6936973769bSHong Zhang MatScalar *aa=a->a,*aval,*ba=b->a,*bval,*ca=c->a,*cval; 694e286af10SHong Zhang Mat_MatMatTransMult *multtrans; 6952128a86cSHong Zhang PetscContainer container; 6966973769bSHong Zhang #if defined(USE_ARRAY) 6976973769bSHong Zhang MatScalar *spdot; 6986973769bSHong Zhang #endif 6995df89d91SHong Zhang 7005df89d91SHong Zhang PetscFunctionBegin; 701e286af10SHong Zhang ierr = PetscObjectQuery((PetscObject)C,"Mat_MatMatTransMult",(PetscObject *)&container);CHKERRQ(ierr); 7022128a86cSHong Zhang if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit"); 7032128a86cSHong Zhang ierr = PetscContainerGetPointer(container,(void **)&multtrans);CHKERRQ(ierr); 7042128a86cSHong Zhang if (multtrans->usecoloring){ 705b9af6bddSHong Zhang MatTransposeColoring matcoloring = multtrans->matcoloring; 706c8db22f6SHong Zhang Mat Bt_dense; 707c8db22f6SHong Zhang PetscInt m,n; 708b2d2b04fSHong Zhang PetscLogDouble t0,tf,etime0=0.0,etime1=0.0,etime2=0.0; 709a0b95ffbSSatish Balay Mat C_dense = multtrans->ABt_den; 710a0b95ffbSSatish Balay 7112128a86cSHong Zhang Bt_dense = multtrans->Bt_den; 712c8db22f6SHong Zhang ierr = MatGetLocalSize(Bt_dense,&m,&n);CHKERRQ(ierr); 713c8db22f6SHong Zhang 714b9af6bddSHong Zhang /* Get Bt_dense by Apply MatTransposeColoring to B */ 7152128a86cSHong Zhang ierr = PetscGetTime(&t0);CHKERRQ(ierr); 716b9af6bddSHong Zhang ierr = MatTransColoringApplySpToDen(matcoloring,B,Bt_dense);CHKERRQ(ierr); 7172128a86cSHong Zhang ierr = PetscGetTime(&tf);CHKERRQ(ierr); 718b2d2b04fSHong Zhang etime0 += tf - t0; 719c8db22f6SHong Zhang 720c8db22f6SHong Zhang /* C_dense = A*Bt_dense */ 7212128a86cSHong Zhang ierr = PetscGetTime(&t0);CHKERRQ(ierr); 7222128a86cSHong Zhang ierr = MatMatMultNumeric_SeqAIJ_SeqDense(A,Bt_dense,C_dense);CHKERRQ(ierr); 7232128a86cSHong Zhang ierr = PetscGetTime(&tf);CHKERRQ(ierr); 7242128a86cSHong Zhang etime2 += tf - t0; 725c8db22f6SHong Zhang 7262128a86cSHong Zhang /* Recover C from C_dense */ 7272128a86cSHong Zhang ierr = PetscGetTime(&t0);CHKERRQ(ierr); 728b9af6bddSHong Zhang ierr = MatTransColoringApplyDenToSp(matcoloring,C_dense,C);CHKERRQ(ierr); 7292128a86cSHong Zhang ierr = PetscGetTime(&tf);CHKERRQ(ierr); 7302128a86cSHong Zhang etime1 += tf - t0; 731f94ccd6cSHong Zhang #if defined(PETSC_USE_INFO) 732f94ccd6cSHong Zhang ierr = PetscInfo4(C,"Num = ColoringApply: %g %g + Mult_sp_dense %g = %g\n",etime0,etime1,etime2,etime0+etime1+etime2); 733f94ccd6cSHong Zhang #endif 734c8db22f6SHong Zhang PetscFunctionReturn(0); 735c8db22f6SHong Zhang } 736c8db22f6SHong Zhang 7376973769bSHong Zhang #if defined(USE_ARRAY) 7386973769bSHong Zhang /* allocate an array for implementing sparse inner-product */ 739e2cac8caSJed Brown ierr = PetscMalloc((A->cmap->n+1)*sizeof(MatScalar),&spdot);CHKERRQ(ierr); 740e2cac8caSJed Brown ierr = PetscMemzero(spdot,(A->cmap->n+1)*sizeof(MatScalar));CHKERRQ(ierr); 7416973769bSHong Zhang #endif 7426973769bSHong Zhang 74381d82fe4SHong Zhang /* clear old values in C */ 74481d82fe4SHong Zhang ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr); 7455df89d91SHong Zhang 7461fa1209cSHong Zhang for (i=0; i<cm; i++) { 74781d82fe4SHong Zhang anzi = ai[i+1] - ai[i]; 74881d82fe4SHong Zhang acol = aj + ai[i]; 7496973769bSHong Zhang aval = aa + ai[i]; 7501fa1209cSHong Zhang cnzi = ci[i+1] - ci[i]; 7511fa1209cSHong Zhang ccol = cj + ci[i]; 7526973769bSHong Zhang cval = ca + ci[i]; 7531fa1209cSHong Zhang for (j=0; j<cnzi; j++){ 75481d82fe4SHong Zhang brow = ccol[j]; 75581d82fe4SHong Zhang bnzj = bi[brow+1] - bi[brow]; 75681d82fe4SHong Zhang bcol = bj + bi[brow]; 7576973769bSHong Zhang bval = ba + bi[brow]; 7586973769bSHong Zhang 75981d82fe4SHong Zhang /* perform sparse inner-product c(i,j)=A[i,:]*B[j,:]^T */ 7606973769bSHong Zhang #if defined(USE_ARRAY) 7616973769bSHong Zhang /* put ba to spdot array */ 7626973769bSHong Zhang for (nextb=0; nextb<bnzj; nextb++) spdot[bcol[nextb]] = bval[nextb]; 7636973769bSHong Zhang /* c(i,j)=A[i,:]*B[j,:]^T */ 7646973769bSHong Zhang for (nexta=0; nexta<anzi; nexta++){ 7656973769bSHong Zhang cval[j] += spdot[acol[nexta]]*aval[nexta]; 7666973769bSHong Zhang } 7676973769bSHong Zhang /* zero spdot array */ 7686973769bSHong Zhang for (nextb=0; nextb<bnzj; nextb++) spdot[bcol[nextb]] = 0.0; 7696973769bSHong Zhang #else 77081d82fe4SHong Zhang nexta = 0; nextb = 0; 77181d82fe4SHong Zhang while (nexta<anzi && nextb<bnzj){ 77281d82fe4SHong Zhang while (acol[nexta] < bcol[nextb] && nexta < anzi) nexta++; 77381d82fe4SHong Zhang if (nexta == anzi) break; 77481d82fe4SHong Zhang while (acol[nexta] > bcol[nextb] && nextb < bnzj) nextb++; 77581d82fe4SHong Zhang if (nextb == bnzj) break; 77681d82fe4SHong Zhang if (acol[nexta] == bcol[nextb]){ 7776973769bSHong Zhang cval[j] += aval[nexta]*bval[nextb]; 77881d82fe4SHong Zhang nexta++; nextb++; 77981d82fe4SHong Zhang flops += 2; 7801fa1209cSHong Zhang } 7811fa1209cSHong Zhang } 7826973769bSHong Zhang #endif 78381d82fe4SHong Zhang } 78481d82fe4SHong Zhang } 78581d82fe4SHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 78681d82fe4SHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 78781d82fe4SHong Zhang ierr = PetscLogFlops(flops);CHKERRQ(ierr); 7886973769bSHong Zhang #if defined(USE_ARRAY) 7896973769bSHong Zhang ierr = PetscFree(spdot); 7906973769bSHong Zhang #endif 7915df89d91SHong Zhang PetscFunctionReturn(0); 7925df89d91SHong Zhang } 7935df89d91SHong Zhang 7945df89d91SHong Zhang #undef __FUNCT__ 79575648e8dSHong Zhang #define __FUNCT__ "MatTransposeMatMult_SeqAIJ_SeqAIJ" 79675648e8dSHong Zhang PetscErrorCode MatTransposeMatMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) { 7975df89d91SHong Zhang PetscErrorCode ierr; 7985df89d91SHong Zhang 7995df89d91SHong Zhang PetscFunctionBegin; 8005df89d91SHong Zhang if (scall == MAT_INITIAL_MATRIX){ 80175648e8dSHong Zhang ierr = MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr); 8025df89d91SHong Zhang } 80375648e8dSHong Zhang ierr = MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr); 8045df89d91SHong Zhang PetscFunctionReturn(0); 8055df89d91SHong Zhang } 8065df89d91SHong Zhang 8075df89d91SHong Zhang #undef __FUNCT__ 80875648e8dSHong Zhang #define __FUNCT__ "MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ" 80975648e8dSHong Zhang PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C) 8105df89d91SHong Zhang { 811bc011b1eSHong Zhang PetscErrorCode ierr; 812bc011b1eSHong Zhang Mat At; 81338baddfdSBarry Smith PetscInt *ati,*atj; 814bc011b1eSHong Zhang 815bc011b1eSHong Zhang PetscFunctionBegin; 816bc011b1eSHong Zhang /* create symbolic At */ 817bc011b1eSHong Zhang ierr = MatGetSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr); 818d0f46423SBarry Smith ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,A->cmap->n,A->rmap->n,ati,atj,PETSC_NULL,&At);CHKERRQ(ierr); 819bc011b1eSHong Zhang 820bc011b1eSHong Zhang /* get symbolic C=At*B */ 821bc011b1eSHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(At,B,fill,C);CHKERRQ(ierr); 822bc011b1eSHong Zhang 823bc011b1eSHong Zhang /* clean up */ 8246bf464f9SBarry Smith ierr = MatDestroy(&At);CHKERRQ(ierr); 825bc011b1eSHong Zhang ierr = MatRestoreSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr); 826bc011b1eSHong Zhang PetscFunctionReturn(0); 827bc011b1eSHong Zhang } 828bc011b1eSHong Zhang 829bc011b1eSHong Zhang #undef __FUNCT__ 83075648e8dSHong Zhang #define __FUNCT__ "MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ" 83175648e8dSHong Zhang PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C) 832bc011b1eSHong Zhang { 833bc011b1eSHong Zhang PetscErrorCode ierr; 8340fbc74f4SHong Zhang Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data; 835d0f46423SBarry Smith PetscInt am=A->rmap->n,anzi,*ai=a->i,*aj=a->j,*bi=b->i,*bj,bnzi,nextb; 836d0f46423SBarry Smith PetscInt cm=C->rmap->n,*ci=c->i,*cj=c->j,crow,*cjj,i,j,k; 837d13dce4bSSatish Balay PetscLogDouble flops=0.0; 8380fbc74f4SHong Zhang MatScalar *aa=a->a,*ba,*ca=c->a,*caj; 839bc011b1eSHong Zhang 840bc011b1eSHong Zhang PetscFunctionBegin; 841bc011b1eSHong Zhang /* clear old values in C */ 842bc011b1eSHong Zhang ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr); 843bc011b1eSHong Zhang 844bc011b1eSHong Zhang /* compute A^T*B using outer product (A^T)[:,i]*B[i,:] */ 845bc011b1eSHong Zhang for (i=0;i<am;i++) { 846bc011b1eSHong Zhang bj = b->j + bi[i]; 847bc011b1eSHong Zhang ba = b->a + bi[i]; 848bc011b1eSHong Zhang bnzi = bi[i+1] - bi[i]; 849bc011b1eSHong Zhang anzi = ai[i+1] - ai[i]; 850bc011b1eSHong Zhang for (j=0; j<anzi; j++) { 851bc011b1eSHong Zhang nextb = 0; 8520fbc74f4SHong Zhang crow = *aj++; 853bc011b1eSHong Zhang cjj = cj + ci[crow]; 854bc011b1eSHong Zhang caj = ca + ci[crow]; 855bc011b1eSHong Zhang /* perform sparse axpy operation. Note cjj includes bj. */ 856bc011b1eSHong Zhang for (k=0; nextb<bnzi; k++) { 8570fbc74f4SHong Zhang if (cjj[k] == *(bj+nextb)) { /* ccol == bcol */ 8580fbc74f4SHong Zhang caj[k] += (*aa)*(*(ba+nextb)); 859bc011b1eSHong Zhang nextb++; 860bc011b1eSHong Zhang } 861bc011b1eSHong Zhang } 862bc011b1eSHong Zhang flops += 2*bnzi; 8630fbc74f4SHong Zhang aa++; 864bc011b1eSHong Zhang } 865bc011b1eSHong Zhang } 866bc011b1eSHong Zhang 867bc011b1eSHong Zhang /* Assemble the final matrix and clean up */ 868bc011b1eSHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 869bc011b1eSHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 870bc011b1eSHong Zhang ierr = PetscLogFlops(flops);CHKERRQ(ierr); 871bc011b1eSHong Zhang PetscFunctionReturn(0); 872bc011b1eSHong Zhang } 8739123193aSHong Zhang 874ec01deb9SMatthew Knepley EXTERN_C_BEGIN 8759123193aSHong Zhang #undef __FUNCT__ 8769123193aSHong Zhang #define __FUNCT__ "MatMatMult_SeqAIJ_SeqDense" 8779123193aSHong Zhang PetscErrorCode MatMatMult_SeqAIJ_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 8789123193aSHong Zhang { 8799123193aSHong Zhang PetscErrorCode ierr; 8809123193aSHong Zhang 8819123193aSHong Zhang PetscFunctionBegin; 8829123193aSHong Zhang if (scall == MAT_INITIAL_MATRIX){ 8839123193aSHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqDense(A,B,fill,C);CHKERRQ(ierr); 8849123193aSHong Zhang } 8859123193aSHong Zhang ierr = MatMatMultNumeric_SeqAIJ_SeqDense(A,B,*C);CHKERRQ(ierr); 8869123193aSHong Zhang PetscFunctionReturn(0); 8879123193aSHong Zhang } 888ec01deb9SMatthew Knepley EXTERN_C_END 889edd81eecSMatthew Knepley 8909123193aSHong Zhang #undef __FUNCT__ 8919123193aSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqDense" 8929123193aSHong Zhang PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C) 8939123193aSHong Zhang { 8949123193aSHong Zhang PetscErrorCode ierr; 8959123193aSHong Zhang 8969123193aSHong Zhang PetscFunctionBegin; 8975a586d82SBarry Smith ierr = MatMatMultSymbolic_SeqDense_SeqDense(A,B,0.0,C);CHKERRQ(ierr); 8988cdbd757SHong Zhang (*C)->ops->matmult = MatMatMult_SeqAIJ_SeqDense; 8999123193aSHong Zhang PetscFunctionReturn(0); 9009123193aSHong Zhang } 9019123193aSHong Zhang 9029123193aSHong Zhang #undef __FUNCT__ 9039123193aSHong Zhang #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqDense" 9049123193aSHong Zhang PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqDense(Mat A,Mat B,Mat C) 9059123193aSHong Zhang { 906f32d5d43SBarry Smith Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data; 9079123193aSHong Zhang PetscErrorCode ierr; 908dd6ea824SBarry Smith PetscScalar *b,*c,r1,r2,r3,r4,*b1,*b2,*b3,*b4; 909dd6ea824SBarry Smith MatScalar *aa; 910d0f46423SBarry Smith PetscInt cm=C->rmap->n, cn=B->cmap->n, bm=B->rmap->n, col, i,j,n,*aj, am = A->rmap->n; 911f32d5d43SBarry Smith PetscInt am2 = 2*am, am3 = 3*am, bm4 = 4*bm,colam; 9129123193aSHong Zhang 9139123193aSHong Zhang PetscFunctionBegin; 914f32d5d43SBarry Smith if (!cm || !cn) PetscFunctionReturn(0); 915e32f2f54SBarry 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); 916e32f2f54SBarry 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); 917e32f2f54SBarry 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); 918f32d5d43SBarry Smith ierr = MatGetArray(B,&b);CHKERRQ(ierr); 919f32d5d43SBarry Smith ierr = MatGetArray(C,&c);CHKERRQ(ierr); 920f32d5d43SBarry Smith b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm; 921f32d5d43SBarry Smith for (col=0; col<cn-4; col += 4){ /* over columns of C */ 922f32d5d43SBarry Smith colam = col*am; 923f32d5d43SBarry Smith for (i=0; i<am; i++) { /* over rows of C in those columns */ 924f32d5d43SBarry Smith r1 = r2 = r3 = r4 = 0.0; 925f32d5d43SBarry Smith n = a->i[i+1] - a->i[i]; 926f32d5d43SBarry Smith aj = a->j + a->i[i]; 927f32d5d43SBarry Smith aa = a->a + a->i[i]; 928f32d5d43SBarry Smith for (j=0; j<n; j++) { 929f32d5d43SBarry Smith r1 += (*aa)*b1[*aj]; 930f32d5d43SBarry Smith r2 += (*aa)*b2[*aj]; 931f32d5d43SBarry Smith r3 += (*aa)*b3[*aj]; 932f32d5d43SBarry Smith r4 += (*aa++)*b4[*aj++]; 9339123193aSHong Zhang } 934f32d5d43SBarry Smith c[colam + i] = r1; 935f32d5d43SBarry Smith c[colam + am + i] = r2; 936f32d5d43SBarry Smith c[colam + am2 + i] = r3; 937f32d5d43SBarry Smith c[colam + am3 + i] = r4; 938f32d5d43SBarry Smith } 939f32d5d43SBarry Smith b1 += bm4; 940f32d5d43SBarry Smith b2 += bm4; 941f32d5d43SBarry Smith b3 += bm4; 942f32d5d43SBarry Smith b4 += bm4; 943f32d5d43SBarry Smith } 944f32d5d43SBarry Smith for (;col<cn; col++){ /* over extra columns of C */ 945f32d5d43SBarry Smith for (i=0; i<am; i++) { /* over rows of C in those columns */ 946f32d5d43SBarry Smith r1 = 0.0; 947f32d5d43SBarry Smith n = a->i[i+1] - a->i[i]; 948f32d5d43SBarry Smith aj = a->j + a->i[i]; 949f32d5d43SBarry Smith aa = a->a + a->i[i]; 950f32d5d43SBarry Smith 951f32d5d43SBarry Smith for (j=0; j<n; j++) { 952f32d5d43SBarry Smith r1 += (*aa++)*b1[*aj++]; 953f32d5d43SBarry Smith } 954f32d5d43SBarry Smith c[col*am + i] = r1; 955f32d5d43SBarry Smith } 956f32d5d43SBarry Smith b1 += bm; 957f32d5d43SBarry Smith } 958dc0b31edSSatish Balay ierr = PetscLogFlops(cn*(2.0*a->nz));CHKERRQ(ierr); 959f32d5d43SBarry Smith ierr = MatRestoreArray(B,&b);CHKERRQ(ierr); 960f32d5d43SBarry Smith ierr = MatRestoreArray(C,&c);CHKERRQ(ierr); 961f32d5d43SBarry Smith ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 962f32d5d43SBarry Smith ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 963f32d5d43SBarry Smith PetscFunctionReturn(0); 964f32d5d43SBarry Smith } 965f32d5d43SBarry Smith 966f32d5d43SBarry Smith /* 9674324174eSBarry Smith Note very similar to MatMult_SeqAIJ(), should generate both codes from same base 968f32d5d43SBarry Smith */ 969f32d5d43SBarry Smith #undef __FUNCT__ 970f32d5d43SBarry Smith #define __FUNCT__ "MatMatMultNumericAdd_SeqAIJ_SeqDense" 971f32d5d43SBarry Smith PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat A,Mat B,Mat C) 972f32d5d43SBarry Smith { 973f32d5d43SBarry Smith Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data; 974f32d5d43SBarry Smith PetscErrorCode ierr; 975dd6ea824SBarry Smith PetscScalar *b,*c,r1,r2,r3,r4,*b1,*b2,*b3,*b4; 976dd6ea824SBarry Smith MatScalar *aa; 977d0f46423SBarry 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; 9784324174eSBarry Smith PetscInt am2 = 2*am, am3 = 3*am, bm4 = 4*bm,colam,*ridx; 979f32d5d43SBarry Smith 980f32d5d43SBarry Smith PetscFunctionBegin; 981f32d5d43SBarry Smith if (!cm || !cn) PetscFunctionReturn(0); 982f32d5d43SBarry Smith ierr = MatGetArray(B,&b);CHKERRQ(ierr); 983f32d5d43SBarry Smith ierr = MatGetArray(C,&c);CHKERRQ(ierr); 984f32d5d43SBarry Smith b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm; 9854324174eSBarry Smith 9864324174eSBarry Smith if (a->compressedrow.use){ /* use compressed row format */ 9874324174eSBarry Smith for (col=0; col<cn-4; col += 4){ /* over columns of C */ 9884324174eSBarry Smith colam = col*am; 9894324174eSBarry Smith arm = a->compressedrow.nrows; 9904324174eSBarry Smith ii = a->compressedrow.i; 9914324174eSBarry Smith ridx = a->compressedrow.rindex; 9924324174eSBarry Smith for (i=0; i<arm; i++) { /* over rows of C in those columns */ 9934324174eSBarry Smith r1 = r2 = r3 = r4 = 0.0; 9944324174eSBarry Smith n = ii[i+1] - ii[i]; 9954324174eSBarry Smith aj = a->j + ii[i]; 9964324174eSBarry Smith aa = a->a + ii[i]; 9974324174eSBarry Smith for (j=0; j<n; j++) { 9984324174eSBarry Smith r1 += (*aa)*b1[*aj]; 9994324174eSBarry Smith r2 += (*aa)*b2[*aj]; 10004324174eSBarry Smith r3 += (*aa)*b3[*aj]; 10014324174eSBarry Smith r4 += (*aa++)*b4[*aj++]; 10024324174eSBarry Smith } 10034324174eSBarry Smith c[colam + ridx[i]] += r1; 10044324174eSBarry Smith c[colam + am + ridx[i]] += r2; 10054324174eSBarry Smith c[colam + am2 + ridx[i]] += r3; 10064324174eSBarry Smith c[colam + am3 + ridx[i]] += r4; 10074324174eSBarry Smith } 10084324174eSBarry Smith b1 += bm4; 10094324174eSBarry Smith b2 += bm4; 10104324174eSBarry Smith b3 += bm4; 10114324174eSBarry Smith b4 += bm4; 10124324174eSBarry Smith } 10134324174eSBarry Smith for (;col<cn; col++){ /* over extra columns of C */ 10144324174eSBarry Smith colam = col*am; 10154324174eSBarry Smith arm = a->compressedrow.nrows; 10164324174eSBarry Smith ii = a->compressedrow.i; 10174324174eSBarry Smith ridx = a->compressedrow.rindex; 10184324174eSBarry Smith for (i=0; i<arm; i++) { /* over rows of C in those columns */ 10194324174eSBarry Smith r1 = 0.0; 10204324174eSBarry Smith n = ii[i+1] - ii[i]; 10214324174eSBarry Smith aj = a->j + ii[i]; 10224324174eSBarry Smith aa = a->a + ii[i]; 10234324174eSBarry Smith 10244324174eSBarry Smith for (j=0; j<n; j++) { 10254324174eSBarry Smith r1 += (*aa++)*b1[*aj++]; 10264324174eSBarry Smith } 10274324174eSBarry Smith c[col*am + ridx[i]] += r1; 10284324174eSBarry Smith } 10294324174eSBarry Smith b1 += bm; 10304324174eSBarry Smith } 10314324174eSBarry Smith } else { 1032f32d5d43SBarry Smith for (col=0; col<cn-4; col += 4){ /* over columns of C */ 1033f32d5d43SBarry Smith colam = col*am; 1034f32d5d43SBarry Smith for (i=0; i<am; i++) { /* over rows of C in those columns */ 1035f32d5d43SBarry Smith r1 = r2 = r3 = r4 = 0.0; 1036f32d5d43SBarry Smith n = a->i[i+1] - a->i[i]; 1037f32d5d43SBarry Smith aj = a->j + a->i[i]; 1038f32d5d43SBarry Smith aa = a->a + a->i[i]; 1039f32d5d43SBarry Smith for (j=0; j<n; j++) { 1040f32d5d43SBarry Smith r1 += (*aa)*b1[*aj]; 1041f32d5d43SBarry Smith r2 += (*aa)*b2[*aj]; 1042f32d5d43SBarry Smith r3 += (*aa)*b3[*aj]; 1043f32d5d43SBarry Smith r4 += (*aa++)*b4[*aj++]; 1044f32d5d43SBarry Smith } 1045f32d5d43SBarry Smith c[colam + i] += r1; 1046f32d5d43SBarry Smith c[colam + am + i] += r2; 1047f32d5d43SBarry Smith c[colam + am2 + i] += r3; 1048f32d5d43SBarry Smith c[colam + am3 + i] += r4; 1049f32d5d43SBarry Smith } 1050f32d5d43SBarry Smith b1 += bm4; 1051f32d5d43SBarry Smith b2 += bm4; 1052f32d5d43SBarry Smith b3 += bm4; 1053f32d5d43SBarry Smith b4 += bm4; 1054f32d5d43SBarry Smith } 1055f32d5d43SBarry Smith for (;col<cn; col++){ /* over extra columns of C */ 1056f32d5d43SBarry Smith for (i=0; i<am; i++) { /* over rows of C in those columns */ 1057f32d5d43SBarry Smith r1 = 0.0; 1058f32d5d43SBarry Smith n = a->i[i+1] - a->i[i]; 1059f32d5d43SBarry Smith aj = a->j + a->i[i]; 1060f32d5d43SBarry Smith aa = a->a + a->i[i]; 1061f32d5d43SBarry Smith 1062f32d5d43SBarry Smith for (j=0; j<n; j++) { 1063f32d5d43SBarry Smith r1 += (*aa++)*b1[*aj++]; 1064f32d5d43SBarry Smith } 1065f32d5d43SBarry Smith c[col*am + i] += r1; 1066f32d5d43SBarry Smith } 1067f32d5d43SBarry Smith b1 += bm; 1068f32d5d43SBarry Smith } 10694324174eSBarry Smith } 1070dc0b31edSSatish Balay ierr = PetscLogFlops(cn*2.0*a->nz);CHKERRQ(ierr); 1071f32d5d43SBarry Smith ierr = MatRestoreArray(B,&b);CHKERRQ(ierr); 1072f32d5d43SBarry Smith ierr = MatRestoreArray(C,&c);CHKERRQ(ierr); 10739123193aSHong Zhang PetscFunctionReturn(0); 10749123193aSHong Zhang } 1075b1683b59SHong Zhang 1076b1683b59SHong Zhang #undef __FUNCT__ 1077b9af6bddSHong Zhang #define __FUNCT__ "MatTransColoringApplySpToDen_SeqAIJ" 1078b9af6bddSHong Zhang PetscErrorCode MatTransColoringApplySpToDen_SeqAIJ(MatTransposeColoring coloring,Mat B,Mat Btdense) 1079c8db22f6SHong Zhang { 1080c8db22f6SHong Zhang PetscErrorCode ierr; 10812f5f1f90SHong Zhang Mat_SeqAIJ *b = (Mat_SeqAIJ*)B->data; 10822f5f1f90SHong Zhang Mat_SeqDense *btdense = (Mat_SeqDense*)Btdense->data; 10832f5f1f90SHong Zhang PetscInt *bi=b->i,*bj=b->j; 10842f5f1f90SHong Zhang PetscInt m=Btdense->rmap->n,n=Btdense->cmap->n,j,k,l,col,anz,*btcol,brow,ncolumns; 10852f5f1f90SHong Zhang MatScalar *btval,*btval_den,*ba=b->a; 10862f5f1f90SHong Zhang PetscInt *columns=coloring->columns,*colorforcol=coloring->colorforcol,ncolors=coloring->ncolors; 1087c8db22f6SHong Zhang 1088c8db22f6SHong Zhang PetscFunctionBegin; 10892f5f1f90SHong Zhang btval_den=btdense->v; 10902f5f1f90SHong Zhang ierr = PetscMemzero(btval_den,(m*n)*sizeof(MatScalar));CHKERRQ(ierr); 10912f5f1f90SHong Zhang for (k=0; k<ncolors; k++) { 10922f5f1f90SHong Zhang ncolumns = coloring->ncolumns[k]; 10932f5f1f90SHong Zhang for (l=0; l<ncolumns; l++) { /* insert a row of B to a column of Btdense */ 1094d2853540SHong Zhang col = *(columns + colorforcol[k] + l); 10952f5f1f90SHong Zhang btcol = bj + bi[col]; 10962f5f1f90SHong Zhang btval = ba + bi[col]; 10972f5f1f90SHong Zhang anz = bi[col+1] - bi[col]; 1098d2853540SHong Zhang for (j=0; j<anz; j++){ 10992f5f1f90SHong Zhang brow = btcol[j]; 11002f5f1f90SHong Zhang btval_den[brow] = btval[j]; 1101c8db22f6SHong Zhang } 1102c8db22f6SHong Zhang } 11032f5f1f90SHong Zhang btval_den += m; 1104c8db22f6SHong Zhang } 1105c8db22f6SHong Zhang PetscFunctionReturn(0); 1106c8db22f6SHong Zhang } 1107c8db22f6SHong Zhang 1108c8db22f6SHong Zhang #undef __FUNCT__ 1109b9af6bddSHong Zhang #define __FUNCT__ "MatTransColoringApplyDenToSp_SeqAIJ" 1110b9af6bddSHong Zhang PetscErrorCode MatTransColoringApplyDenToSp_SeqAIJ(MatTransposeColoring matcoloring,Mat Cden,Mat Csp) 11118972f759SHong Zhang { 11128972f759SHong Zhang PetscErrorCode ierr; 1113b2d2b04fSHong Zhang Mat_SeqAIJ *csp = (Mat_SeqAIJ*)Csp->data; 11142f5f1f90SHong Zhang PetscInt k,l,*row,*idx,m,ncolors=matcoloring->ncolors,nrows; 1115b2d2b04fSHong Zhang PetscScalar *ca_den,*cp_den,*ca=csp->a; 11162f5f1f90SHong Zhang PetscInt *rows=matcoloring->rows,*spidx=matcoloring->columnsforspidx,*colorforrow=matcoloring->colorforrow; 11178972f759SHong Zhang 11188972f759SHong Zhang PetscFunctionBegin; 11198972f759SHong Zhang ierr = MatGetLocalSize(Csp,&m,PETSC_NULL);CHKERRQ(ierr); 1120b2d2b04fSHong Zhang ierr = MatGetArray(Cden,&ca_den);CHKERRQ(ierr); 1121b2d2b04fSHong Zhang cp_den = ca_den; 11222f5f1f90SHong Zhang for (k=0; k<ncolors; k++) { 11232f5f1f90SHong Zhang nrows = matcoloring->nrows[k]; 11242f5f1f90SHong Zhang row = rows + colorforrow[k]; 11252f5f1f90SHong Zhang idx = spidx + colorforrow[k]; 11262f5f1f90SHong Zhang for (l=0; l<nrows; l++){ 11272f5f1f90SHong Zhang ca[idx[l]] = cp_den[row[l]]; 1128b2d2b04fSHong Zhang } 1129b2d2b04fSHong Zhang cp_den += m; 1130b2d2b04fSHong Zhang } 1131b2d2b04fSHong Zhang ierr = MatRestoreArray(Cden,&ca_den);CHKERRQ(ierr); 11328972f759SHong Zhang PetscFunctionReturn(0); 11338972f759SHong Zhang } 11348972f759SHong Zhang 1135*e99be685SHong Zhang /* 1136*e99be685SHong Zhang MatGetColumnIJ_SeqAIJ_Color() and MatRestoreColumnIJ_SeqAIJ_Color() are customized from 1137*e99be685SHong Zhang MatGetColumnIJ_SeqAIJ() and MatRestoreColumnIJ_SeqAIJ() by adding an output 1138*e99be685SHong Zhang spidx[], index of a->j, to be used for setting 'columnsforspidx' in MatTransposeColoringCreate_SeqAIJ(). 1139*e99be685SHong Zhang */ 1140*e99be685SHong Zhang #undef __FUNCT__ 1141*e99be685SHong Zhang #define __FUNCT__ "MatGetColumnIJ_SeqAIJ_Color" 1142*e99be685SHong Zhang PetscErrorCode MatGetColumnIJ_SeqAIJ_Color(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool inodecompressed,PetscInt *nn,PetscInt *ia[],PetscInt *ja[],PetscInt *spidx[],PetscBool *done) 1143*e99be685SHong Zhang { 1144*e99be685SHong Zhang Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 1145*e99be685SHong Zhang PetscErrorCode ierr; 1146*e99be685SHong Zhang PetscInt i,*collengths,*cia,*cja,n = A->cmap->n,m = A->rmap->n; 1147*e99be685SHong Zhang PetscInt nz = a->i[m],row,*jj,mr,col; 1148*e99be685SHong Zhang PetscInt *cspidx; 1149*e99be685SHong Zhang 1150*e99be685SHong Zhang PetscFunctionBegin; 1151*e99be685SHong Zhang *nn = n; 1152*e99be685SHong Zhang if (!ia) PetscFunctionReturn(0); 1153*e99be685SHong Zhang if (symmetric) { 1154*e99be685SHong Zhang SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"MatGetColumnIJ_SeqAIJ_Color() not supported for the case symmetric"); 1155*e99be685SHong Zhang ierr = MatToSymmetricIJ_SeqAIJ(A->rmap->n,a->i,a->j,0,oshift,ia,ja);CHKERRQ(ierr); 1156*e99be685SHong Zhang } else { 1157*e99be685SHong Zhang ierr = PetscMalloc((n+1)*sizeof(PetscInt),&collengths);CHKERRQ(ierr); 1158*e99be685SHong Zhang ierr = PetscMemzero(collengths,n*sizeof(PetscInt));CHKERRQ(ierr); 1159*e99be685SHong Zhang ierr = PetscMalloc((n+1)*sizeof(PetscInt),&cia);CHKERRQ(ierr); 1160*e99be685SHong Zhang ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&cja);CHKERRQ(ierr); 1161*e99be685SHong Zhang ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&cspidx);CHKERRQ(ierr); 1162*e99be685SHong Zhang jj = a->j; 1163*e99be685SHong Zhang for (i=0; i<nz; i++) { 1164*e99be685SHong Zhang collengths[jj[i]]++; 1165*e99be685SHong Zhang } 1166*e99be685SHong Zhang cia[0] = oshift; 1167*e99be685SHong Zhang for (i=0; i<n; i++) { 1168*e99be685SHong Zhang cia[i+1] = cia[i] + collengths[i]; 1169*e99be685SHong Zhang } 1170*e99be685SHong Zhang ierr = PetscMemzero(collengths,n*sizeof(PetscInt));CHKERRQ(ierr); 1171*e99be685SHong Zhang jj = a->j; 1172*e99be685SHong Zhang for (row=0; row<m; row++) { 1173*e99be685SHong Zhang mr = a->i[row+1] - a->i[row]; 1174*e99be685SHong Zhang for (i=0; i<mr; i++) { 1175*e99be685SHong Zhang col = *jj++; 1176*e99be685SHong Zhang cspidx[cia[col] + collengths[col] - oshift] = a->i[row] + i; /* index of a->j */ 1177*e99be685SHong Zhang cja[cia[col] + collengths[col]++ - oshift] = row + oshift; 1178*e99be685SHong Zhang } 1179*e99be685SHong Zhang } 1180*e99be685SHong Zhang ierr = PetscFree(collengths);CHKERRQ(ierr); 1181*e99be685SHong Zhang *ia = cia; *ja = cja; 1182*e99be685SHong Zhang *spidx = cspidx; 1183*e99be685SHong Zhang } 1184*e99be685SHong Zhang PetscFunctionReturn(0); 1185*e99be685SHong Zhang } 1186*e99be685SHong Zhang 1187*e99be685SHong Zhang #undef __FUNCT__ 1188*e99be685SHong Zhang #define __FUNCT__ "MatRestoreColumnIJ_SeqAIJ_Color" 1189*e99be685SHong Zhang PetscErrorCode MatRestoreColumnIJ_SeqAIJ_Color(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool inodecompressed,PetscInt *n,PetscInt *ia[],PetscInt *ja[],PetscInt *spidx[],PetscBool *done) 1190*e99be685SHong Zhang { 1191*e99be685SHong Zhang PetscErrorCode ierr; 1192*e99be685SHong Zhang 1193*e99be685SHong Zhang PetscFunctionBegin; 1194*e99be685SHong Zhang if (!ia) PetscFunctionReturn(0); 1195*e99be685SHong Zhang 1196*e99be685SHong Zhang ierr = PetscFree(*ia);CHKERRQ(ierr); 1197*e99be685SHong Zhang ierr = PetscFree(*ja);CHKERRQ(ierr); 1198*e99be685SHong Zhang ierr = PetscFree(*spidx);CHKERRQ(ierr); 1199*e99be685SHong Zhang PetscFunctionReturn(0); 1200*e99be685SHong Zhang } 1201*e99be685SHong Zhang 12028972f759SHong Zhang #undef __FUNCT__ 1203b9af6bddSHong Zhang #define __FUNCT__ "MatTransposeColoringCreate_SeqAIJ" 1204b9af6bddSHong Zhang PetscErrorCode MatTransposeColoringCreate_SeqAIJ(Mat mat,ISColoring iscoloring,MatTransposeColoring c) 1205b1683b59SHong Zhang { 1206b1683b59SHong Zhang PetscErrorCode ierr; 1207d2853540SHong Zhang PetscInt i,n,nrows,N,j,k,m,*row_idx,*ci,*cj,ncols,col,cm; 1208b1683b59SHong Zhang const PetscInt *is; 1209b28d80bdSHong Zhang PetscInt nis = iscoloring->n,*rowhit,bs = 1; 1210b1683b59SHong Zhang IS *isa; 1211b1683b59SHong Zhang PetscBool done; 1212b1683b59SHong Zhang PetscBool flg1,flg2; 1213b28d80bdSHong Zhang Mat_SeqAIJ *csp = (Mat_SeqAIJ*)mat->data; 1214*e99be685SHong Zhang PetscInt *colorforrow,*rows,*rows_i,*columnsforspidx,*columnsforspidx_i,*idxhit,*spidx; 1215d2853540SHong Zhang PetscInt *colorforcol,*columns,*columns_i; 1216b1683b59SHong Zhang 1217b1683b59SHong Zhang PetscFunctionBegin; 1218b1683b59SHong Zhang ierr = ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);CHKERRQ(ierr); 1219*e99be685SHong Zhang 1220b1683b59SHong Zhang /* this is ugly way to get blocksize but cannot call MatGetBlockSize() because AIJ can have bs > 1 */ 1221b1683b59SHong Zhang ierr = PetscTypeCompare((PetscObject)mat,MATSEQBAIJ,&flg1);CHKERRQ(ierr); 1222b1683b59SHong Zhang ierr = PetscTypeCompare((PetscObject)mat,MATMPIBAIJ,&flg2);CHKERRQ(ierr); 1223b1683b59SHong Zhang if (flg1 || flg2) { 1224b1683b59SHong Zhang ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr); 1225b1683b59SHong Zhang } 1226b1683b59SHong Zhang 1227b1683b59SHong Zhang N = mat->cmap->N/bs; 1228b1683b59SHong Zhang c->M = mat->rmap->N/bs; /* set total rows, columns and local rows */ 1229b1683b59SHong Zhang c->N = mat->cmap->N/bs; 1230b1683b59SHong Zhang c->m = mat->rmap->N/bs; 1231b1683b59SHong Zhang c->rstart = 0; 1232b1683b59SHong Zhang 1233b1683b59SHong Zhang c->ncolors = nis; 1234b1683b59SHong Zhang ierr = PetscMalloc(nis*sizeof(PetscInt),&c->ncolumns);CHKERRQ(ierr); 1235b28d80bdSHong Zhang ierr = PetscMalloc(nis*sizeof(PetscInt),&c->nrows);CHKERRQ(ierr); 1236d2853540SHong Zhang ierr = PetscMalloc2(csp->nz+1,PetscInt,&rows,csp->nz+1,PetscInt,&columnsforspidx);CHKERRQ(ierr); 1237d2853540SHong Zhang ierr = PetscMalloc((nis+1)*sizeof(PetscInt),&colorforrow);CHKERRQ(ierr); 1238d2853540SHong Zhang colorforrow[0] = 0; 1239d2853540SHong Zhang rows_i = rows; 1240d2853540SHong Zhang columnsforspidx_i = columnsforspidx; 1241b1683b59SHong Zhang 1242d2853540SHong Zhang ierr = PetscMalloc((nis+1)*sizeof(PetscInt),&colorforcol);CHKERRQ(ierr); 1243d2853540SHong Zhang ierr = PetscMalloc((N+1)*sizeof(PetscInt),&columns);CHKERRQ(ierr); 1244d2853540SHong Zhang colorforcol[0] = 0; 1245d2853540SHong Zhang columns_i = columns; 1246d2853540SHong Zhang 1247*e99be685SHong Zhang ierr = MatGetColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,&done);CHKERRQ(ierr); /* column-wise storage of mat */ 1248b1683b59SHong Zhang if (!done) SETERRQ1(((PetscObject)mat)->comm,PETSC_ERR_SUP,"MatGetColumnIJ() not supported for matrix type %s",((PetscObject)mat)->type_name); 1249b1683b59SHong Zhang 1250b28d80bdSHong Zhang cm = c->m; 1251b28d80bdSHong Zhang ierr = PetscMalloc((cm+1)*sizeof(PetscInt),&rowhit);CHKERRQ(ierr); 1252*e99be685SHong Zhang ierr = PetscMalloc((cm+1)*sizeof(PetscInt),&idxhit);CHKERRQ(ierr); 1253b1683b59SHong Zhang for (i=0; i<nis; i++) { 1254b1683b59SHong Zhang ierr = ISGetLocalSize(isa[i],&n);CHKERRQ(ierr); 1255b1683b59SHong Zhang ierr = ISGetIndices(isa[i],&is);CHKERRQ(ierr); 1256b1683b59SHong Zhang c->ncolumns[i] = n; 1257b1683b59SHong Zhang if (n) { 1258d2853540SHong Zhang ierr = PetscMemcpy(columns_i,is,n*sizeof(PetscInt));CHKERRQ(ierr); 1259b1683b59SHong Zhang } 1260d2853540SHong Zhang colorforcol[i+1] = colorforcol[i] + n; 1261d2853540SHong Zhang columns_i += n; 1262b1683b59SHong Zhang 1263b1683b59SHong Zhang /* fast, crude version requires O(N*N) work */ 1264e8354b3eSHong Zhang ierr = PetscMemzero(rowhit,cm*sizeof(PetscInt));CHKERRQ(ierr); 1265*e99be685SHong Zhang 1266b1683b59SHong Zhang /* loop over columns*/ 1267b1683b59SHong Zhang for (j=0; j<n; j++) { 1268b1683b59SHong Zhang col = is[j]; 1269d2853540SHong Zhang row_idx = cj + ci[col]; 1270b1683b59SHong Zhang m = ci[col+1] - ci[col]; 1271b1683b59SHong Zhang /* loop over columns marking them in rowhit */ 1272b1683b59SHong Zhang for (k=0; k<m; k++) { 1273*e99be685SHong Zhang idxhit[*row_idx] = spidx[ci[col] + k]; 1274d2853540SHong Zhang rowhit[*row_idx++] = col + 1; 1275b1683b59SHong Zhang } 1276b1683b59SHong Zhang } 1277b1683b59SHong Zhang /* count the number of hits */ 1278b1683b59SHong Zhang nrows = 0; 1279e8354b3eSHong Zhang for (j=0; j<cm; j++) { 1280b1683b59SHong Zhang if (rowhit[j]) nrows++; 1281b1683b59SHong Zhang } 1282b1683b59SHong Zhang c->nrows[i] = nrows; 1283d2853540SHong Zhang colorforrow[i+1] = colorforrow[i] + nrows; 1284d2853540SHong Zhang 1285b1683b59SHong Zhang nrows = 0; 1286e8354b3eSHong Zhang for (j=0; j<cm; j++) { 1287b1683b59SHong Zhang if (rowhit[j]) { 1288d2853540SHong Zhang rows_i[nrows] = j; 1289*e99be685SHong Zhang columnsforspidx_i[nrows] = idxhit[j]; 1290b1683b59SHong Zhang nrows++; 1291b1683b59SHong Zhang } 1292b1683b59SHong Zhang } 1293b1683b59SHong Zhang ierr = ISRestoreIndices(isa[i],&is);CHKERRQ(ierr); 1294d2853540SHong Zhang rows_i += nrows; columnsforspidx_i += nrows; 1295b1683b59SHong Zhang } 1296*e99be685SHong Zhang ierr = MatRestoreColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,&done);CHKERRQ(ierr); 1297b28d80bdSHong Zhang ierr = PetscFree(rowhit);CHKERRQ(ierr); 1298b1683b59SHong Zhang ierr = ISColoringRestoreIS(iscoloring,&isa);CHKERRQ(ierr); 1299d2853540SHong Zhang #if defined(PETSC_USE_DEBUG) 1300d2853540SHong Zhang if (csp->nz != colorforrow[nis]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"csp->nz %d != colorforrow[nis] %d",csp->nz,colorforrow[nis]); 1301d2853540SHong Zhang #endif 1302b28d80bdSHong Zhang 1303d2853540SHong Zhang c->colorforrow = colorforrow; 1304d2853540SHong Zhang c->rows = rows; 1305d2853540SHong Zhang c->columnsforspidx = columnsforspidx; 1306d2853540SHong Zhang c->colorforcol = colorforcol; 1307d2853540SHong Zhang c->columns = columns; 1308f94ccd6cSHong Zhang ierr = PetscFree(idxhit);CHKERRQ(ierr); 1309b1683b59SHong Zhang PetscFunctionReturn(0); 1310b1683b59SHong Zhang } 1311