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; 203c50cad2SHong Zhang PetscBool scalable=PETSC_FALSE,scalable_fast=PETSC_FALSE; 215c66b693SKris Buschelman 225c66b693SKris Buschelman PetscFunctionBegin; 2326be0446SHong Zhang if (scall == MAT_INITIAL_MATRIX){ 24152983bfSHong Zhang ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr); 25152983bfSHong Zhang ierr = PetscOptionsBool("-matmatmult_scalable","Use a scalable but slower C=A*B","",scalable,&scalable,PETSC_NULL);CHKERRQ(ierr); 263c50cad2SHong Zhang ierr = PetscOptionsBool("-matmatmult_scalable_fast","Use a scalable but slower C=A*B","",scalable_fast,&scalable_fast,PETSC_NULL);CHKERRQ(ierr); 27d8bbc50fSBarry Smith ierr = PetscOptionsEnd();CHKERRQ(ierr); 28d8bbc50fSBarry Smith ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr); 293c50cad2SHong Zhang if (scalable_fast){ 303c50cad2SHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast(A,B,fill,C);CHKERRQ(ierr); 31aacf854cSBarry Smith } else if (scalable){ 32aacf854cSBarry Smith ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(A,B,fill,C);CHKERRQ(ierr); 3325023028SHong Zhang } else { 3426be0446SHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr); 3525023028SHong Zhang } 36598bc09dSHong Zhang ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr); 3726be0446SHong Zhang } 38598bc09dSHong Zhang ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr); 3901e47db4SHong Zhang ierr = (*(*C)->ops->matmultnumeric)(A,B,*C);CHKERRQ(ierr); 40598bc09dSHong Zhang ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr); 415c66b693SKris Buschelman PetscFunctionReturn(0); 425c66b693SKris Buschelman } 43ec01deb9SMatthew Knepley EXTERN_C_END 441c24bd37SHong Zhang 45e9e4536cSHong Zhang #undef __FUNCT__ 46b561aa9dSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ" 47b561aa9dSHong Zhang PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C) 48b561aa9dSHong Zhang { 49b561aa9dSHong Zhang PetscErrorCode ierr; 50b561aa9dSHong Zhang Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c; 51971236abSHong Zhang PetscInt *ai=a->i,*bi=b->i,*ci,*cj; 52*a2f3521dSMark F. Adams PetscInt am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N,rbs=A->rmap->bs,cbs=B->cmap->bs; 53b561aa9dSHong Zhang PetscReal afill; 54fb908031SHong Zhang PetscInt i,j,anzi,brow,bnzj,cnzi,*bj,*aj,nlnk_max,*lnk,ndouble=0; 55fb908031SHong Zhang PetscBT lnkbt; 56fb908031SHong Zhang PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 57b561aa9dSHong Zhang 58b561aa9dSHong Zhang PetscFunctionBegin; 59bd958071SHong Zhang /* Get ci and cj */ 60bd958071SHong Zhang /*---------------*/ 61fb908031SHong Zhang /* Allocate ci array, arrays for fill computation and */ 62fb908031SHong Zhang /* free space for accumulating nonzero column info */ 63fb908031SHong Zhang ierr = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr); 64fb908031SHong Zhang ci[0] = 0; 65fb908031SHong Zhang 66fb908031SHong Zhang /* create and initialize a linked list */ 67fb908031SHong Zhang nlnk_max = a->rmax*b->rmax; 68fb908031SHong Zhang if (!nlnk_max || nlnk_max > bn) nlnk_max = bn; 69fb908031SHong Zhang ierr = PetscLLCondensedCreate(nlnk_max,bn,&lnk,&lnkbt);CHKERRQ(ierr); 70fb908031SHong Zhang 71fb908031SHong Zhang /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */ 72fb908031SHong Zhang ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr); 73fb908031SHong Zhang current_space = free_space; 74fb908031SHong Zhang 75fb908031SHong Zhang /* Determine ci and cj */ 76fb908031SHong Zhang for (i=0; i<am; i++) { 77fb908031SHong Zhang anzi = ai[i+1] - ai[i]; 78fb908031SHong Zhang aj = a->j + ai[i]; 79fb908031SHong Zhang for (j=0; j<anzi; j++){ 80fb908031SHong Zhang brow = aj[j]; 81fb908031SHong Zhang bnzj = bi[brow+1] - bi[brow]; 82fb908031SHong Zhang bj = b->j + bi[brow]; 83fb908031SHong Zhang /* add non-zero cols of B into the sorted linked list lnk */ 84fb908031SHong Zhang ierr = PetscLLCondensedAddSorted(bnzj,bj,lnk,lnkbt);CHKERRQ(ierr); 85fb908031SHong Zhang } 86fb908031SHong Zhang cnzi = lnk[0]; 87fb908031SHong Zhang 88fb908031SHong Zhang /* If free space is not available, make more free space */ 89fb908031SHong Zhang /* Double the amount of total space in the list */ 90fb908031SHong Zhang if (current_space->local_remaining<cnzi) { 91fb908031SHong Zhang ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 92fb908031SHong Zhang ndouble++; 93fb908031SHong Zhang } 94fb908031SHong Zhang 95fb908031SHong Zhang /* Copy data into free space, then initialize lnk */ 96fb908031SHong Zhang ierr = PetscLLCondensedClean(bn,cnzi,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 97fb908031SHong Zhang current_space->array += cnzi; 98fb908031SHong Zhang current_space->local_used += cnzi; 99fb908031SHong Zhang current_space->local_remaining -= cnzi; 100fb908031SHong Zhang ci[i+1] = ci[i] + cnzi; 101fb908031SHong Zhang } 102fb908031SHong Zhang 103fb908031SHong Zhang /* Column indices are in the list of free space */ 104fb908031SHong Zhang /* Allocate space for cj, initialize cj, and */ 105fb908031SHong Zhang /* destroy list of free space and other temporary array(s) */ 106fb908031SHong Zhang ierr = PetscMalloc((ci[am]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr); 107fb908031SHong Zhang ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr); 108fb908031SHong Zhang ierr = PetscLLCondensedDestroy(lnk,lnkbt);CHKERRQ(ierr); 109b561aa9dSHong Zhang 11026be0446SHong Zhang /* put together the new symbolic matrix */ 111aa1aec99SHong Zhang ierr = MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,am,bn,ci,cj,PETSC_NULL,C);CHKERRQ(ierr); 112*a2f3521dSMark F. Adams (*C)->rmap->bs = A->rmap->bs; 113*a2f3521dSMark F. Adams (*C)->cmap->bs = B->cmap->bs; 11458c24d83SHong Zhang 11558c24d83SHong Zhang /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 11658c24d83SHong Zhang /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */ 11758c24d83SHong Zhang c = (Mat_SeqAIJ *)((*C)->data); 118aa1aec99SHong Zhang c->free_a = PETSC_FALSE; 119e6b907acSBarry Smith c->free_ij = PETSC_TRUE; 12058c24d83SHong Zhang c->nonew = 0; 121002e8affSHong Zhang (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ; /* fast, needs non-scalable O(bn) array 'abdense' */ 1220b7e3e3dSHong Zhang 1237212da7cSHong Zhang /* set MatInfo */ 1247212da7cSHong Zhang afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5; 125f2b054eeSHong Zhang if (afill < 1.0) afill = 1.0; 1267212da7cSHong Zhang c->maxnz = ci[am]; 1277212da7cSHong Zhang c->nz = ci[am]; 128fb908031SHong Zhang (*C)->info.mallocs = ndouble; 1297212da7cSHong Zhang (*C)->info.fill_ratio_given = fill; 1307212da7cSHong Zhang (*C)->info.fill_ratio_needed = afill; 1317212da7cSHong Zhang 1327212da7cSHong Zhang #if defined(PETSC_USE_INFO) 1337212da7cSHong Zhang if (ci[am]) { 134fb908031SHong Zhang ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",ndouble,fill,afill);CHKERRQ(ierr); 135f2b054eeSHong Zhang ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr); 136f2b054eeSHong Zhang } else { 137f2b054eeSHong Zhang ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr); 138be0fcf8dSHong Zhang } 139f2b054eeSHong Zhang #endif 14058c24d83SHong Zhang PetscFunctionReturn(0); 14158c24d83SHong Zhang } 142d50806bdSBarry Smith 14326be0446SHong Zhang #undef __FUNCT__ 14426be0446SHong Zhang #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqAIJ" 145dfbe8321SBarry Smith PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C) 146d50806bdSBarry Smith { 147dfbe8321SBarry Smith PetscErrorCode ierr; 148d13dce4bSSatish Balay PetscLogDouble flops=0.0; 149d50806bdSBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 150d50806bdSBarry Smith Mat_SeqAIJ *b = (Mat_SeqAIJ *)B->data; 151d50806bdSBarry Smith Mat_SeqAIJ *c = (Mat_SeqAIJ *)C->data; 15238baddfdSBarry Smith PetscInt *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j; 153c58ca2e3SHong Zhang PetscInt am=A->rmap->n,cm=C->rmap->n; 154519eb980SHong Zhang PetscInt i,j,k,anzi,bnzi,cnzi,brow; 155aa1aec99SHong Zhang PetscScalar *aa=a->a,*ba=b->a,*baj,*ca,valtmp; 156aa1aec99SHong Zhang PetscScalar *ab_dense; 157d50806bdSBarry Smith 158d50806bdSBarry Smith PetscFunctionBegin; 159aa1aec99SHong Zhang /* printf("MatMatMultNumeric_SeqAIJ_SeqAIJ...ca %p\n",c->a); */ 160aa1aec99SHong Zhang if (!c->a){ /* first call of MatMatMultNumeric_SeqAIJ_SeqAIJ, allocate ca and matmult_abdense */ 161aa1aec99SHong Zhang ierr = PetscMalloc((ci[cm]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr); 162aa1aec99SHong Zhang c->a = ca; 163aa1aec99SHong Zhang c->free_a = PETSC_TRUE; 164aa1aec99SHong Zhang 165aa1aec99SHong Zhang ierr = PetscMalloc(B->cmap->N*sizeof(PetscScalar),&ab_dense);CHKERRQ(ierr); 166aa1aec99SHong Zhang ierr = PetscMemzero(ab_dense,B->cmap->N*sizeof(PetscScalar));CHKERRQ(ierr); 167aa1aec99SHong Zhang c->matmult_abdense = ab_dense; 168aa1aec99SHong Zhang } else { 169aa1aec99SHong Zhang ca = c->a; 170aa1aec99SHong Zhang ab_dense = c->matmult_abdense; 171aa1aec99SHong Zhang } 172aa1aec99SHong Zhang 173c124e916SHong Zhang /* clean old values in C */ 174c124e916SHong Zhang ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr); 175d50806bdSBarry Smith /* Traverse A row-wise. */ 176d50806bdSBarry Smith /* Build the ith row in C by summing over nonzero columns in A, */ 177d50806bdSBarry Smith /* the rows of B corresponding to nonzeros of A. */ 178d50806bdSBarry Smith for (i=0; i<am; i++) { 179d50806bdSBarry Smith anzi = ai[i+1] - ai[i]; 180d50806bdSBarry Smith for (j=0; j<anzi; j++) { 181519eb980SHong Zhang brow = aj[j]; 182d50806bdSBarry Smith bnzi = bi[brow+1] - bi[brow]; 183d50806bdSBarry Smith bjj = bj + bi[brow]; 184d50806bdSBarry Smith baj = ba + bi[brow]; 185519eb980SHong Zhang /* perform dense axpy */ 18636ec6d2dSHong Zhang valtmp = aa[j]; 187519eb980SHong Zhang for (k=0; k<bnzi; k++) { 18836ec6d2dSHong Zhang ab_dense[bjj[k]] += valtmp*baj[k]; 189519eb980SHong Zhang } 190519eb980SHong Zhang flops += 2*bnzi; 191519eb980SHong Zhang } 192c58ca2e3SHong Zhang aj += anzi; aa += anzi; 193c58ca2e3SHong Zhang 194c58ca2e3SHong Zhang cnzi = ci[i+1] - ci[i]; 195519eb980SHong Zhang for (k=0; k<cnzi; k++) { 196519eb980SHong Zhang ca[k] += ab_dense[cj[k]]; 197519eb980SHong Zhang ab_dense[cj[k]] = 0.0; /* zero ab_dense */ 198519eb980SHong Zhang } 199519eb980SHong Zhang flops += cnzi; 200519eb980SHong Zhang cj += cnzi; ca += cnzi; 201519eb980SHong Zhang } 202c58ca2e3SHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 203c58ca2e3SHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 204c58ca2e3SHong Zhang ierr = PetscLogFlops(flops);CHKERRQ(ierr); 205c58ca2e3SHong Zhang PetscFunctionReturn(0); 206c58ca2e3SHong Zhang } 207c58ca2e3SHong Zhang 208c58ca2e3SHong Zhang #undef __FUNCT__ 20925023028SHong Zhang #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable" 21025023028SHong Zhang PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(Mat A,Mat B,Mat C) 211c58ca2e3SHong Zhang { 212c58ca2e3SHong Zhang PetscErrorCode ierr; 213c58ca2e3SHong Zhang PetscLogDouble flops=0.0; 214c58ca2e3SHong Zhang Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 215c58ca2e3SHong Zhang Mat_SeqAIJ *b = (Mat_SeqAIJ *)B->data; 216c58ca2e3SHong Zhang Mat_SeqAIJ *c = (Mat_SeqAIJ *)C->data; 217c58ca2e3SHong Zhang PetscInt *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j; 218c58ca2e3SHong Zhang PetscInt am=A->rmap->N,cm=C->rmap->N; 219c58ca2e3SHong Zhang PetscInt i,j,k,anzi,bnzi,cnzi,brow; 22036ec6d2dSHong Zhang PetscScalar *aa=a->a,*ba=b->a,*baj,*ca=c->a,valtmp; 221c58ca2e3SHong Zhang PetscInt nextb; 222c58ca2e3SHong Zhang 223c58ca2e3SHong Zhang PetscFunctionBegin; 224e9e4536cSHong Zhang #if defined(DEBUG_MATMATMULT) 225971236abSHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF,"MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable...\n");CHKERRQ(ierr); 226e9e4536cSHong Zhang #endif 227c58ca2e3SHong Zhang /* clean old values in C */ 228c58ca2e3SHong Zhang ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr); 229c58ca2e3SHong Zhang /* Traverse A row-wise. */ 230c58ca2e3SHong Zhang /* Build the ith row in C by summing over nonzero columns in A, */ 231c58ca2e3SHong Zhang /* the rows of B corresponding to nonzeros of A. */ 232519eb980SHong Zhang for (i=0;i<am;i++) { 233519eb980SHong Zhang anzi = ai[i+1] - ai[i]; 234519eb980SHong Zhang cnzi = ci[i+1] - ci[i]; 235519eb980SHong Zhang for (j=0;j<anzi;j++) { 236519eb980SHong Zhang brow = aj[j]; 237519eb980SHong Zhang bnzi = bi[brow+1] - bi[brow]; 238519eb980SHong Zhang bjj = bj + bi[brow]; 239519eb980SHong Zhang baj = ba + bi[brow]; 240519eb980SHong Zhang /* perform sparse axpy */ 24136ec6d2dSHong Zhang valtmp = aa[j]; 242c124e916SHong Zhang nextb = 0; 243c124e916SHong Zhang for (k=0; nextb<bnzi; k++) { 244c124e916SHong Zhang if (cj[k] == bjj[nextb]){ /* ccol == bcol */ 24536ec6d2dSHong Zhang ca[k] += valtmp*baj[nextb++]; 246c124e916SHong Zhang } 247d50806bdSBarry Smith } 248d50806bdSBarry Smith flops += 2*bnzi; 249d50806bdSBarry Smith } 250519eb980SHong Zhang aj += anzi; aa += anzi; 251519eb980SHong Zhang cj += cnzi; ca += cnzi; 252d50806bdSBarry Smith } 253c58ca2e3SHong Zhang 254716bacf3SKris Buschelman ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 255716bacf3SKris Buschelman ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 256d50806bdSBarry Smith ierr = PetscLogFlops(flops);CHKERRQ(ierr); 257d50806bdSBarry Smith PetscFunctionReturn(0); 258d50806bdSBarry Smith } 259bc011b1eSHong Zhang 260e9e4536cSHong Zhang #undef __FUNCT__ 2613c50cad2SHong Zhang #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast" 2623c50cad2SHong Zhang PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast(Mat A,Mat B,PetscReal fill,Mat *C) 26325296bd5SBarry Smith { 26425296bd5SBarry Smith PetscErrorCode ierr; 26525296bd5SBarry Smith Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c; 26625296bd5SBarry Smith PetscInt *ai=a->i,*bi=b->i,*ci,*cj; 2673c50cad2SHong Zhang PetscInt am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N; 26825296bd5SBarry Smith MatScalar *ca; 26925296bd5SBarry Smith PetscReal afill; 2703c50cad2SHong Zhang PetscInt i,j,anzi,brow,bnzj,cnzi,*bj,*aj,nlnk_max,*lnk,ndouble=0; 2713c50cad2SHong Zhang PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 27225296bd5SBarry Smith 27325296bd5SBarry Smith PetscFunctionBegin; 2743c50cad2SHong Zhang /* Get ci and cj - same as MatMatMultSymbolic_SeqAIJ_SeqAIJ except using PetscLLxxx_fast() */ 2753c50cad2SHong Zhang /*-----------------------------------------------------------------------------------------*/ 2763c50cad2SHong Zhang /* Allocate arrays for fill computation and free space for accumulating nonzero column */ 2773c50cad2SHong Zhang ierr = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr); 2783c50cad2SHong Zhang ci[0] = 0; 2793c50cad2SHong Zhang 2803c50cad2SHong Zhang /* create and initialize a linked list */ 2813c50cad2SHong Zhang nlnk_max = a->rmax*b->rmax; 2823c50cad2SHong Zhang if (!nlnk_max || nlnk_max > bn) nlnk_max = bn; /* in case rmax is not defined for A or B */ 2833c50cad2SHong Zhang ierr = PetscLLCondensedCreate_fast(nlnk_max,&lnk);CHKERRQ(ierr); 2843c50cad2SHong Zhang 2853c50cad2SHong Zhang /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */ 2863c50cad2SHong Zhang ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr); 2873c50cad2SHong Zhang current_space = free_space; 2883c50cad2SHong Zhang 2893c50cad2SHong Zhang /* Determine ci and cj */ 2903c50cad2SHong Zhang for (i=0; i<am; i++) { 2913c50cad2SHong Zhang anzi = ai[i+1] - ai[i]; 2923c50cad2SHong Zhang aj = a->j + ai[i]; 2933c50cad2SHong Zhang for (j=0; j<anzi; j++){ 2943c50cad2SHong Zhang brow = aj[j]; 2953c50cad2SHong Zhang bnzj = bi[brow+1] - bi[brow]; 2963c50cad2SHong Zhang bj = b->j + bi[brow]; 2973c50cad2SHong Zhang /* add non-zero cols of B into the sorted linked list lnk */ 2983c50cad2SHong Zhang ierr = PetscLLCondensedAddSorted_fast(bnzj,bj,lnk);CHKERRQ(ierr); 2993c50cad2SHong Zhang } 3003c50cad2SHong Zhang cnzi = lnk[1]; 3013c50cad2SHong Zhang 3023c50cad2SHong Zhang /* If free space is not available, make more free space */ 3033c50cad2SHong Zhang /* Double the amount of total space in the list */ 3043c50cad2SHong Zhang if (current_space->local_remaining<cnzi) { 3053c50cad2SHong Zhang ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 3063c50cad2SHong Zhang ndouble++; 3073c50cad2SHong Zhang } 3083c50cad2SHong Zhang 3093c50cad2SHong Zhang /* Copy data into free space, then initialize lnk */ 3103c50cad2SHong Zhang ierr = PetscLLCondensedClean_fast(cnzi,current_space->array,lnk);CHKERRQ(ierr); 3113c50cad2SHong Zhang current_space->array += cnzi; 3123c50cad2SHong Zhang current_space->local_used += cnzi; 3133c50cad2SHong Zhang current_space->local_remaining -= cnzi; 3143c50cad2SHong Zhang ci[i+1] = ci[i] + cnzi; 3153c50cad2SHong Zhang } 3163c50cad2SHong Zhang 3173c50cad2SHong Zhang /* Column indices are in the list of free space */ 3183c50cad2SHong Zhang /* Allocate space for cj, initialize cj, and */ 3193c50cad2SHong Zhang /* destroy list of free space and other temporary array(s) */ 3203c50cad2SHong Zhang ierr = PetscMalloc((ci[am]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr); 3213c50cad2SHong Zhang ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr); 3223c50cad2SHong Zhang ierr = PetscLLCondensedDestroy_fast(lnk);CHKERRQ(ierr); 32325296bd5SBarry Smith 32425296bd5SBarry Smith /* Allocate space for ca */ 32525296bd5SBarry Smith ierr = PetscMalloc((ci[am]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr); 32625296bd5SBarry Smith ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr); 32725296bd5SBarry Smith 32825296bd5SBarry Smith /* put together the new symbolic matrix */ 32925296bd5SBarry Smith ierr = MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,am,bn,ci,cj,ca,C);CHKERRQ(ierr); 330*a2f3521dSMark F. Adams (*C)->rmap->bs = A->rmap->bs; 331*a2f3521dSMark F. Adams (*C)->cmap->bs = B->cmap->bs; 33225296bd5SBarry Smith 33325296bd5SBarry Smith /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 33425296bd5SBarry Smith /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */ 33525296bd5SBarry Smith c = (Mat_SeqAIJ *)((*C)->data); 33625296bd5SBarry Smith c->free_a = PETSC_TRUE; 33725296bd5SBarry Smith c->free_ij = PETSC_TRUE; 33825296bd5SBarry Smith c->nonew = 0; 33925296bd5SBarry Smith (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable; /* slower, less memory */ 34025296bd5SBarry Smith 34125296bd5SBarry Smith /* set MatInfo */ 34225296bd5SBarry Smith afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5; 34325296bd5SBarry Smith if (afill < 1.0) afill = 1.0; 34425296bd5SBarry Smith c->maxnz = ci[am]; 34525296bd5SBarry Smith c->nz = ci[am]; 3463c50cad2SHong Zhang (*C)->info.mallocs = ndouble; 34725296bd5SBarry Smith (*C)->info.fill_ratio_given = fill; 34825296bd5SBarry Smith (*C)->info.fill_ratio_needed = afill; 34925296bd5SBarry Smith 35025296bd5SBarry Smith #if defined(PETSC_USE_INFO) 35125296bd5SBarry Smith if (ci[am]) { 3523c50cad2SHong Zhang ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",ndouble,fill,afill);CHKERRQ(ierr); 35325296bd5SBarry Smith ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr); 35425296bd5SBarry Smith } else { 35525296bd5SBarry Smith ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr); 35625296bd5SBarry Smith } 35725296bd5SBarry Smith #endif 35825296bd5SBarry Smith PetscFunctionReturn(0); 35925296bd5SBarry Smith } 36025296bd5SBarry Smith 36125296bd5SBarry Smith 36225296bd5SBarry Smith #undef __FUNCT__ 36325023028SHong Zhang #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable" 36425023028SHong Zhang PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(Mat A,Mat B,PetscReal fill,Mat *C) 365e9e4536cSHong Zhang { 366e9e4536cSHong Zhang PetscErrorCode ierr; 367e9e4536cSHong Zhang Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c; 368bf9555e6SHong Zhang PetscInt *ai=a->i,*bi=b->i,*ci,*cj; 36925c41797SHong Zhang PetscInt am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N; 370e9e4536cSHong Zhang MatScalar *ca; 371e9e4536cSHong Zhang PetscReal afill; 372bd958071SHong Zhang PetscInt i,j,anzi,brow,bnzj,cnzi,*bj,*aj,nlnk_max,*lnk,ndouble=0; 373bd958071SHong Zhang PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 374e9e4536cSHong Zhang 375e9e4536cSHong Zhang PetscFunctionBegin; 376bd958071SHong Zhang /* Get ci and cj - same as MatMatMultSymbolic_SeqAIJ_SeqAIJ except using PetscLLxxx_Scalalbe() */ 377bd958071SHong Zhang /*---------------------------------------------------------------------------------------------*/ 378bd958071SHong Zhang /* Allocate arrays for fill computation and free space for accumulating nonzero column */ 379bd958071SHong Zhang ierr = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr); 380bd958071SHong Zhang ci[0] = 0; 381bd958071SHong Zhang 382bd958071SHong Zhang /* create and initialize a linked list */ 383bd958071SHong Zhang nlnk_max = a->rmax*b->rmax; 384bd958071SHong Zhang if (!nlnk_max || nlnk_max > bn) nlnk_max = bn; /* in case rmax is not defined for A or B */ 385bd958071SHong Zhang ierr = PetscLLCondensedCreate_Scalable(nlnk_max,&lnk);CHKERRQ(ierr); 386bd958071SHong Zhang 387bd958071SHong Zhang /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */ 388bd958071SHong Zhang ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr); 389bd958071SHong Zhang current_space = free_space; 390bd958071SHong Zhang 391bd958071SHong Zhang /* Determine ci and cj */ 392bd958071SHong Zhang for (i=0; i<am; i++) { 393bd958071SHong Zhang anzi = ai[i+1] - ai[i]; 394bd958071SHong Zhang aj = a->j + ai[i]; 395bd958071SHong Zhang for (j=0; j<anzi; j++){ 396bd958071SHong Zhang brow = aj[j]; 397bd958071SHong Zhang bnzj = bi[brow+1] - bi[brow]; 398bd958071SHong Zhang bj = b->j + bi[brow]; 399bd958071SHong Zhang /* add non-zero cols of B into the sorted linked list lnk */ 400bd958071SHong Zhang ierr = PetscLLCondensedAddSorted_Scalable(bnzj,bj,lnk);CHKERRQ(ierr); 401bd958071SHong Zhang } 402bd958071SHong Zhang cnzi = lnk[0]; 403bd958071SHong Zhang 404bd958071SHong Zhang /* If free space is not available, make more free space */ 405bd958071SHong Zhang /* Double the amount of total space in the list */ 406bd958071SHong Zhang if (current_space->local_remaining<cnzi) { 407bd958071SHong Zhang ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 408bd958071SHong Zhang ndouble++; 409bd958071SHong Zhang } 410bd958071SHong Zhang 411bd958071SHong Zhang /* Copy data into free space, then initialize lnk */ 412bd958071SHong Zhang ierr = PetscLLCondensedClean_Scalable(cnzi,current_space->array,lnk);CHKERRQ(ierr); 413bd958071SHong Zhang current_space->array += cnzi; 414bd958071SHong Zhang current_space->local_used += cnzi; 415bd958071SHong Zhang current_space->local_remaining -= cnzi; 416bd958071SHong Zhang ci[i+1] = ci[i] + cnzi; 417bd958071SHong Zhang } 418bd958071SHong Zhang 419bd958071SHong Zhang /* Column indices are in the list of free space */ 420bd958071SHong Zhang /* Allocate space for cj, initialize cj, and */ 421bd958071SHong Zhang /* destroy list of free space and other temporary array(s) */ 422bd958071SHong Zhang ierr = PetscMalloc((ci[am]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr); 423bd958071SHong Zhang ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr); 424bd958071SHong Zhang ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr); 425e9e4536cSHong Zhang 426e9e4536cSHong Zhang /* Allocate space for ca */ 427bd958071SHong Zhang /*-----------------------*/ 428e9e4536cSHong Zhang ierr = PetscMalloc((ci[am]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr); 429e9e4536cSHong Zhang ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr); 430e9e4536cSHong Zhang 431e9e4536cSHong Zhang /* put together the new symbolic matrix */ 432e9e4536cSHong Zhang ierr = MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,am,bn,ci,cj,ca,C);CHKERRQ(ierr); 433*a2f3521dSMark F. Adams (*C)->rmap->bs = A->rmap->bs; 434*a2f3521dSMark F. Adams (*C)->cmap->bs = B->cmap->bs; 435e9e4536cSHong Zhang 436e9e4536cSHong Zhang /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 437e9e4536cSHong Zhang /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */ 438e9e4536cSHong Zhang c = (Mat_SeqAIJ *)((*C)->data); 439e9e4536cSHong Zhang c->free_a = PETSC_TRUE; 440e9e4536cSHong Zhang c->free_ij = PETSC_TRUE; 441e9e4536cSHong Zhang c->nonew = 0; 44225023028SHong Zhang (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable; /* slower, less memory */ 443e9e4536cSHong Zhang 444e9e4536cSHong Zhang /* set MatInfo */ 445e9e4536cSHong Zhang afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5; 446e9e4536cSHong Zhang if (afill < 1.0) afill = 1.0; 447e9e4536cSHong Zhang c->maxnz = ci[am]; 448e9e4536cSHong Zhang c->nz = ci[am]; 449bd958071SHong Zhang (*C)->info.mallocs = ndouble; 450e9e4536cSHong Zhang (*C)->info.fill_ratio_given = fill; 451e9e4536cSHong Zhang (*C)->info.fill_ratio_needed = afill; 452e9e4536cSHong Zhang 453e9e4536cSHong Zhang #if defined(PETSC_USE_INFO) 454e9e4536cSHong Zhang if (ci[am]) { 455bd958071SHong Zhang ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",ndouble,fill,afill);CHKERRQ(ierr); 456e9e4536cSHong Zhang ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr); 457e9e4536cSHong Zhang } else { 458e9e4536cSHong Zhang ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr); 459e9e4536cSHong Zhang } 460e9e4536cSHong Zhang #endif 461e9e4536cSHong Zhang PetscFunctionReturn(0); 462e9e4536cSHong Zhang } 463e9e4536cSHong Zhang 464e9e4536cSHong Zhang 465d2853540SHong Zhang /* This routine is not used. Should be removed! */ 466bc011b1eSHong Zhang #undef __FUNCT__ 4676fc122caSHong Zhang #define __FUNCT__ "MatMatTransposeMult_SeqAIJ_SeqAIJ" 4686fc122caSHong Zhang PetscErrorCode MatMatTransposeMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 4695df89d91SHong Zhang { 470bc011b1eSHong Zhang PetscErrorCode ierr; 471bc011b1eSHong Zhang 472bc011b1eSHong Zhang PetscFunctionBegin; 473bc011b1eSHong Zhang if (scall == MAT_INITIAL_MATRIX){ 4746fc122caSHong Zhang ierr = MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr); 475bc011b1eSHong Zhang } 4766fc122caSHong Zhang ierr = MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr); 477bc011b1eSHong Zhang PetscFunctionReturn(0); 478bc011b1eSHong Zhang } 479bc011b1eSHong Zhang 480bc011b1eSHong Zhang #undef __FUNCT__ 481e286af10SHong Zhang #define __FUNCT__ "PetscContainerDestroy_Mat_MatMatTransMult" 482e286af10SHong Zhang PetscErrorCode PetscContainerDestroy_Mat_MatMatTransMult(void *ptr) 4832128a86cSHong Zhang { 4842128a86cSHong Zhang PetscErrorCode ierr; 485e286af10SHong Zhang Mat_MatMatTransMult *multtrans=(Mat_MatMatTransMult*)ptr; 4862128a86cSHong Zhang 4872128a86cSHong Zhang PetscFunctionBegin; 488b9af6bddSHong Zhang ierr = MatTransposeColoringDestroy(&multtrans->matcoloring);CHKERRQ(ierr); 4892128a86cSHong Zhang ierr = MatDestroy(&multtrans->Bt_den);CHKERRQ(ierr); 4902128a86cSHong Zhang ierr = MatDestroy(&multtrans->ABt_den);CHKERRQ(ierr); 4912128a86cSHong Zhang ierr = PetscFree(multtrans);CHKERRQ(ierr); 4922128a86cSHong Zhang PetscFunctionReturn(0); 4932128a86cSHong Zhang } 4942128a86cSHong Zhang 4952128a86cSHong Zhang #undef __FUNCT__ 4962128a86cSHong Zhang #define __FUNCT__ "MatDestroy_SeqAIJ_MatMatMultTrans" 4972128a86cSHong Zhang PetscErrorCode MatDestroy_SeqAIJ_MatMatMultTrans(Mat A) 4982128a86cSHong Zhang { 4992128a86cSHong Zhang PetscErrorCode ierr; 5002128a86cSHong Zhang PetscContainer container; 501e286af10SHong Zhang Mat_MatMatTransMult *multtrans=PETSC_NULL; 5022128a86cSHong Zhang 5032128a86cSHong Zhang PetscFunctionBegin; 504e286af10SHong Zhang ierr = PetscObjectQuery((PetscObject)A,"Mat_MatMatTransMult",(PetscObject *)&container);CHKERRQ(ierr); 5052128a86cSHong Zhang if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit"); 5062128a86cSHong Zhang ierr = PetscContainerGetPointer(container,(void **)&multtrans);CHKERRQ(ierr); 5072128a86cSHong Zhang A->ops->destroy = multtrans->destroy; 5082128a86cSHong Zhang if (A->ops->destroy) { 5092128a86cSHong Zhang ierr = (*A->ops->destroy)(A);CHKERRQ(ierr); 5102128a86cSHong Zhang } 511e286af10SHong Zhang ierr = PetscObjectCompose((PetscObject)A,"Mat_MatMatTransMult",0);CHKERRQ(ierr); 5122128a86cSHong Zhang PetscFunctionReturn(0); 5132128a86cSHong Zhang } 5142128a86cSHong Zhang 5152128a86cSHong Zhang #undef __FUNCT__ 5166fc122caSHong Zhang #define __FUNCT__ "MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ" 5176fc122caSHong Zhang PetscErrorCode MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C) 518bc011b1eSHong Zhang { 5195df89d91SHong Zhang PetscErrorCode ierr; 52081d82fe4SHong Zhang Mat Bt; 52181d82fe4SHong Zhang PetscInt *bti,*btj; 522e286af10SHong Zhang Mat_MatMatTransMult *multtrans; 5232128a86cSHong Zhang PetscContainer container; 524d2853540SHong Zhang PetscLogDouble t0,tf,etime2=0.0; 525d2853540SHong Zhang 52681d82fe4SHong Zhang PetscFunctionBegin; 527d2853540SHong Zhang ierr = PetscGetTime(&t0);CHKERRQ(ierr); 52881d82fe4SHong Zhang /* create symbolic Bt */ 52981d82fe4SHong Zhang ierr = MatGetSymbolicTranspose_SeqAIJ(B,&bti,&btj);CHKERRQ(ierr); 53081d82fe4SHong Zhang ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,B->cmap->n,B->rmap->n,bti,btj,PETSC_NULL,&Bt);CHKERRQ(ierr); 531*a2f3521dSMark F. Adams (*C)->rmap->bs = A->cmap->bs; 532*a2f3521dSMark F. Adams (*C)->cmap->bs = B->cmap->bs; 53381d82fe4SHong Zhang 53481d82fe4SHong Zhang /* get symbolic C=A*Bt */ 53581d82fe4SHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,Bt,fill,C);CHKERRQ(ierr); 53681d82fe4SHong Zhang 5372128a86cSHong Zhang /* create a supporting struct for reuse intermidiate dense matrices with matcoloring */ 538e286af10SHong Zhang ierr = PetscNew(Mat_MatMatTransMult,&multtrans);CHKERRQ(ierr); 5392128a86cSHong Zhang 5402128a86cSHong Zhang /* attach the supporting struct to C */ 5412128a86cSHong Zhang ierr = PetscContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr); 5422128a86cSHong Zhang ierr = PetscContainerSetPointer(container,multtrans);CHKERRQ(ierr); 543e286af10SHong Zhang ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_Mat_MatMatTransMult);CHKERRQ(ierr); 544e286af10SHong Zhang ierr = PetscObjectCompose((PetscObject)(*C),"Mat_MatMatTransMult",(PetscObject)container);CHKERRQ(ierr); 5452128a86cSHong Zhang ierr = PetscContainerDestroy(&container);CHKERRQ(ierr); 5462128a86cSHong Zhang 5472128a86cSHong Zhang multtrans->usecoloring = PETSC_FALSE; 5482128a86cSHong Zhang multtrans->destroy = (*C)->ops->destroy; 5492128a86cSHong Zhang (*C)->ops->destroy = MatDestroy_SeqAIJ_MatMatMultTrans; 5502128a86cSHong Zhang 551d2853540SHong Zhang ierr = PetscGetTime(&tf);CHKERRQ(ierr); 552d2853540SHong Zhang etime2 += tf - t0; 553d2853540SHong Zhang 554f400d4cbSHong Zhang ierr = PetscOptionsGetBool(PETSC_NULL,"-matmattransmult_color",&multtrans->usecoloring,PETSC_NULL);CHKERRQ(ierr); 5552128a86cSHong Zhang if (multtrans->usecoloring){ 556b9af6bddSHong Zhang /* Create MatTransposeColoring from symbolic C=A*B^T */ 557b9af6bddSHong Zhang MatTransposeColoring matcoloring; 5582128a86cSHong Zhang ISColoring iscoloring; 5592128a86cSHong Zhang Mat Bt_dense,C_dense; 560d2853540SHong Zhang PetscLogDouble etime0=0.0,etime01=0.0,etime1=0.0; 561e8354b3eSHong Zhang 562e8354b3eSHong Zhang ierr = PetscGetTime(&t0);CHKERRQ(ierr); 563502de53fSHong Zhang ierr = MatGetColoring(*C,MATCOLORINGLF,&iscoloring);CHKERRQ(ierr); 564d2853540SHong Zhang ierr = PetscGetTime(&tf);CHKERRQ(ierr); 565d2853540SHong Zhang etime0 += tf - t0; 566d2853540SHong Zhang 567d2853540SHong Zhang ierr = PetscGetTime(&t0);CHKERRQ(ierr); 568b9af6bddSHong Zhang ierr = MatTransposeColoringCreate(*C,iscoloring,&matcoloring);CHKERRQ(ierr); 5692128a86cSHong Zhang multtrans->matcoloring = matcoloring; 5702128a86cSHong Zhang ierr = ISColoringDestroy(&iscoloring);CHKERRQ(ierr); 571e8354b3eSHong Zhang ierr = PetscGetTime(&tf);CHKERRQ(ierr); 572d2853540SHong Zhang etime01 += tf - t0; 5732128a86cSHong Zhang 574e8354b3eSHong Zhang ierr = PetscGetTime(&t0);CHKERRQ(ierr); 5752128a86cSHong Zhang /* Create Bt_dense and C_dense = A*Bt_dense */ 5762128a86cSHong Zhang ierr = MatCreate(PETSC_COMM_SELF,&Bt_dense);CHKERRQ(ierr); 5772128a86cSHong Zhang ierr = MatSetSizes(Bt_dense,A->cmap->n,matcoloring->ncolors,A->cmap->n,matcoloring->ncolors);CHKERRQ(ierr); 5782128a86cSHong Zhang ierr = MatSetType(Bt_dense,MATSEQDENSE);CHKERRQ(ierr); 5792128a86cSHong Zhang ierr = MatSeqDenseSetPreallocation(Bt_dense,PETSC_NULL);CHKERRQ(ierr); 5802128a86cSHong Zhang Bt_dense->assembled = PETSC_TRUE; 5812128a86cSHong Zhang multtrans->Bt_den = Bt_dense; 5822128a86cSHong Zhang 5832128a86cSHong Zhang ierr = MatCreate(PETSC_COMM_SELF,&C_dense);CHKERRQ(ierr); 5842128a86cSHong Zhang ierr = MatSetSizes(C_dense,A->rmap->n,matcoloring->ncolors,A->rmap->n,matcoloring->ncolors);CHKERRQ(ierr); 5852128a86cSHong Zhang ierr = MatSetType(C_dense,MATSEQDENSE);CHKERRQ(ierr); 5862128a86cSHong Zhang ierr = MatSeqDenseSetPreallocation(C_dense,PETSC_NULL);CHKERRQ(ierr); 5872128a86cSHong Zhang Bt_dense->assembled = PETSC_TRUE; 5882128a86cSHong Zhang multtrans->ABt_den = C_dense; 589e8354b3eSHong Zhang ierr = PetscGetTime(&tf);CHKERRQ(ierr); 590e8354b3eSHong Zhang etime1 += tf - t0; 591f94ccd6cSHong Zhang 592f94ccd6cSHong Zhang #if defined(PETSC_USE_INFO) 5931ce71dffSSatish Balay { 594f94ccd6cSHong Zhang Mat_SeqAIJ *c=(Mat_SeqAIJ*)(*C)->data; 595f94ccd6cSHong 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)); 596f94ccd6cSHong Zhang ierr = PetscInfo5(*C,"Sym = GetColor %g + ColorCreate %g + MatDenseCreate %g + non-colorSym %g = %g\n",etime0,etime01,etime1,etime2,etime0+etime01+etime1+etime2); 5971ce71dffSSatish Balay } 598f94ccd6cSHong Zhang #endif 5992128a86cSHong Zhang } 600e99be685SHong Zhang /* clean up */ 601e99be685SHong Zhang ierr = MatDestroy(&Bt);CHKERRQ(ierr); 602e99be685SHong Zhang ierr = MatRestoreSymbolicTranspose_SeqAIJ(B,&bti,&btj);CHKERRQ(ierr); 6032128a86cSHong Zhang 604f94ccd6cSHong Zhang 605f94ccd6cSHong Zhang 60681d82fe4SHong Zhang #if defined(INEFFICIENT_ALGORITHM) 60781d82fe4SHong Zhang /* The algorithm below computes am*bm sparse inner-product - inefficient! It will be deleted later. */ 6081fa1209cSHong Zhang PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 6091fa1209cSHong Zhang Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c; 6101fa1209cSHong Zhang PetscInt *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*ci,*cj,*acol,*bcol; 6111fa1209cSHong Zhang PetscInt am=A->rmap->N,bm=B->rmap->N; 6121fa1209cSHong Zhang PetscInt i,j,anzi,bnzj,cnzi,nlnk,*lnk,nspacedouble=0,ka,kb,index[1]; 6131fa1209cSHong Zhang MatScalar *ca; 6141fa1209cSHong Zhang PetscBT lnkbt; 6151fa1209cSHong Zhang PetscReal afill; 6165df89d91SHong Zhang 6171fa1209cSHong Zhang /* Allocate row pointer array ci */ 6181fa1209cSHong Zhang ierr = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr); 6191fa1209cSHong Zhang ci[0] = 0; 6201fa1209cSHong Zhang 6211fa1209cSHong Zhang /* Create and initialize a linked list for C columns */ 6221fa1209cSHong Zhang nlnk = bm+1; 6231fa1209cSHong Zhang ierr = PetscLLCreate(bm,bm,nlnk,lnk,lnkbt);CHKERRQ(ierr); 6241fa1209cSHong Zhang 6251fa1209cSHong Zhang /* Initial FreeSpace with size fill*(nnz(A)+nnz(B)) */ 6261fa1209cSHong Zhang ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr); 6271fa1209cSHong Zhang current_space = free_space; 6281fa1209cSHong Zhang 6291fa1209cSHong Zhang /* Determine symbolic info for each row of the product A*B^T: */ 6301fa1209cSHong Zhang for (i=0; i<am; i++) { 6311fa1209cSHong Zhang anzi = ai[i+1] - ai[i]; 6321fa1209cSHong Zhang cnzi = 0; 6331fa1209cSHong Zhang acol = aj + ai[i]; 6341fa1209cSHong Zhang for (j=0; j<bm; j++){ 6351fa1209cSHong Zhang bnzj = bi[j+1] - bi[j]; 6361fa1209cSHong Zhang bcol= bj + bi[j]; 63781d82fe4SHong Zhang /* sparse inner-product c(i,j)=A[i,:]*B[j,:]^T */ 6381fa1209cSHong Zhang ka = 0; kb = 0; 6391fa1209cSHong Zhang while (ka < anzi && kb < bnzj){ 64081d82fe4SHong Zhang while (acol[ka] < bcol[kb] && ka < anzi) ka++; 64181d82fe4SHong Zhang if (ka == anzi) break; 64281d82fe4SHong Zhang while (acol[ka] > bcol[kb] && kb < bnzj) kb++; 64381d82fe4SHong Zhang if (kb == bnzj) break; 6441fa1209cSHong Zhang if (acol[ka] == bcol[kb]){ /* add nonzero c(i,j) to lnk */ 6451fa1209cSHong Zhang index[0] = j; 6461fa1209cSHong Zhang ierr = PetscLLAdd(1,index,bm,nlnk,lnk,lnkbt);CHKERRQ(ierr); 6471fa1209cSHong Zhang cnzi++; 6481fa1209cSHong Zhang break; 6491fa1209cSHong Zhang } 6501fa1209cSHong Zhang } 6511fa1209cSHong Zhang } 6521fa1209cSHong Zhang 6531fa1209cSHong Zhang /* If free space is not available, make more free space */ 6541fa1209cSHong Zhang /* Double the amount of total space in the list */ 6551fa1209cSHong Zhang if (current_space->local_remaining<cnzi) { 6561fa1209cSHong Zhang ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 6571fa1209cSHong Zhang nspacedouble++; 6581fa1209cSHong Zhang } 6591fa1209cSHong Zhang 6601fa1209cSHong Zhang /* Copy data into free space, then initialize lnk */ 6611fa1209cSHong Zhang ierr = PetscLLClean(bm,bm,cnzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 6621fa1209cSHong Zhang current_space->array += cnzi; 6631fa1209cSHong Zhang current_space->local_used += cnzi; 6641fa1209cSHong Zhang current_space->local_remaining -= cnzi; 6651fa1209cSHong Zhang 6661fa1209cSHong Zhang ci[i+1] = ci[i] + cnzi; 6671fa1209cSHong Zhang } 6681fa1209cSHong Zhang 6691fa1209cSHong Zhang 6701fa1209cSHong Zhang /* Column indices are in the list of free space. 6711fa1209cSHong Zhang Allocate array cj, copy column indices to cj, and destroy list of free space */ 6721fa1209cSHong Zhang ierr = PetscMalloc((ci[am]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr); 6731fa1209cSHong Zhang ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr); 6741fa1209cSHong Zhang ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 6751fa1209cSHong Zhang 6761fa1209cSHong Zhang /* Allocate space for ca */ 6771fa1209cSHong Zhang ierr = PetscMalloc((ci[am]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr); 6781fa1209cSHong Zhang ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr); 6791fa1209cSHong Zhang 6801fa1209cSHong Zhang /* put together the new symbolic matrix */ 6811fa1209cSHong Zhang ierr = MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,am,bm,ci,cj,ca,C);CHKERRQ(ierr); 682*a2f3521dSMark F. Adams (*C)->rmap->bs = A->cmap->bs; 683*a2f3521dSMark F. Adams (*C)->cmap->bs = B->cmap->bs; 6841fa1209cSHong Zhang 6851fa1209cSHong Zhang /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 6861fa1209cSHong Zhang /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */ 6871fa1209cSHong Zhang c = (Mat_SeqAIJ *)((*C)->data); 6881fa1209cSHong Zhang c->free_a = PETSC_TRUE; 6891fa1209cSHong Zhang c->free_ij = PETSC_TRUE; 6901fa1209cSHong Zhang c->nonew = 0; 6911fa1209cSHong Zhang 6921fa1209cSHong Zhang /* set MatInfo */ 6931fa1209cSHong Zhang afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5; 6941fa1209cSHong Zhang if (afill < 1.0) afill = 1.0; 6951fa1209cSHong Zhang c->maxnz = ci[am]; 6961fa1209cSHong Zhang c->nz = ci[am]; 6971fa1209cSHong Zhang (*C)->info.mallocs = nspacedouble; 6981fa1209cSHong Zhang (*C)->info.fill_ratio_given = fill; 6991fa1209cSHong Zhang (*C)->info.fill_ratio_needed = afill; 7001fa1209cSHong Zhang 7011fa1209cSHong Zhang #if defined(PETSC_USE_INFO) 7021fa1209cSHong Zhang if (ci[am]) { 7031fa1209cSHong Zhang ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr); 7046fc122caSHong Zhang ierr = PetscInfo1((*C),"Use MatMatTransposeMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr); 7051fa1209cSHong Zhang } else { 7061fa1209cSHong Zhang ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr); 7071fa1209cSHong Zhang } 7081fa1209cSHong Zhang #endif 70981d82fe4SHong Zhang #endif 7105df89d91SHong Zhang PetscFunctionReturn(0); 7115df89d91SHong Zhang } 7125df89d91SHong Zhang 7136973769bSHong Zhang /* #define USE_ARRAY - for sparse dot product. Slower than !USE_ARRAY */ 7145df89d91SHong Zhang #undef __FUNCT__ 7156fc122caSHong Zhang #define __FUNCT__ "MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ" 7166fc122caSHong Zhang PetscErrorCode MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C) 7175df89d91SHong Zhang { 7185df89d91SHong Zhang PetscErrorCode ierr; 7195df89d91SHong Zhang Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data; 720e2cac8caSJed Brown PetscInt *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,anzi,bnzj,nexta,nextb,*acol,*bcol,brow; 72181d82fe4SHong Zhang PetscInt cm=C->rmap->n,*ci=c->i,*cj=c->j,i,j,cnzi,*ccol; 7225df89d91SHong Zhang PetscLogDouble flops=0.0; 723aa1aec99SHong Zhang MatScalar *aa=a->a,*aval,*ba=b->a,*bval,*ca,*cval; 724e286af10SHong Zhang Mat_MatMatTransMult *multtrans; 7252128a86cSHong Zhang PetscContainer container; 7266973769bSHong Zhang #if defined(USE_ARRAY) 7276973769bSHong Zhang MatScalar *spdot; 7286973769bSHong Zhang #endif 7295df89d91SHong Zhang 7305df89d91SHong Zhang PetscFunctionBegin; 731e286af10SHong Zhang ierr = PetscObjectQuery((PetscObject)C,"Mat_MatMatTransMult",(PetscObject *)&container);CHKERRQ(ierr); 7322128a86cSHong Zhang if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit"); 7332128a86cSHong Zhang ierr = PetscContainerGetPointer(container,(void **)&multtrans);CHKERRQ(ierr); 7342128a86cSHong Zhang if (multtrans->usecoloring){ 735b9af6bddSHong Zhang MatTransposeColoring matcoloring = multtrans->matcoloring; 736c8db22f6SHong Zhang Mat Bt_dense; 737c8db22f6SHong Zhang PetscInt m,n; 738b2d2b04fSHong Zhang PetscLogDouble t0,tf,etime0=0.0,etime1=0.0,etime2=0.0; 739a0b95ffbSSatish Balay Mat C_dense = multtrans->ABt_den; 740a0b95ffbSSatish Balay 7412128a86cSHong Zhang Bt_dense = multtrans->Bt_den; 742c8db22f6SHong Zhang ierr = MatGetLocalSize(Bt_dense,&m,&n);CHKERRQ(ierr); 743c8db22f6SHong Zhang 744b9af6bddSHong Zhang /* Get Bt_dense by Apply MatTransposeColoring to B */ 7452128a86cSHong Zhang ierr = PetscGetTime(&t0);CHKERRQ(ierr); 746b9af6bddSHong Zhang ierr = MatTransColoringApplySpToDen(matcoloring,B,Bt_dense);CHKERRQ(ierr); 7472128a86cSHong Zhang ierr = PetscGetTime(&tf);CHKERRQ(ierr); 748b2d2b04fSHong Zhang etime0 += tf - t0; 749c8db22f6SHong Zhang 750c8db22f6SHong Zhang /* C_dense = A*Bt_dense */ 7512128a86cSHong Zhang ierr = PetscGetTime(&t0);CHKERRQ(ierr); 7522128a86cSHong Zhang ierr = MatMatMultNumeric_SeqAIJ_SeqDense(A,Bt_dense,C_dense);CHKERRQ(ierr); 7532128a86cSHong Zhang ierr = PetscGetTime(&tf);CHKERRQ(ierr); 7542128a86cSHong Zhang etime2 += tf - t0; 755c8db22f6SHong Zhang 7562128a86cSHong Zhang /* Recover C from C_dense */ 7572128a86cSHong Zhang ierr = PetscGetTime(&t0);CHKERRQ(ierr); 758b9af6bddSHong Zhang ierr = MatTransColoringApplyDenToSp(matcoloring,C_dense,C);CHKERRQ(ierr); 7592128a86cSHong Zhang ierr = PetscGetTime(&tf);CHKERRQ(ierr); 7602128a86cSHong Zhang etime1 += tf - t0; 761f94ccd6cSHong Zhang #if defined(PETSC_USE_INFO) 762f94ccd6cSHong Zhang ierr = PetscInfo4(C,"Num = ColoringApply: %g %g + Mult_sp_dense %g = %g\n",etime0,etime1,etime2,etime0+etime1+etime2); 763f94ccd6cSHong Zhang #endif 764c8db22f6SHong Zhang PetscFunctionReturn(0); 765c8db22f6SHong Zhang } 766c8db22f6SHong Zhang 7676973769bSHong Zhang #if defined(USE_ARRAY) 7686973769bSHong Zhang /* allocate an array for implementing sparse inner-product */ 769e2cac8caSJed Brown ierr = PetscMalloc((A->cmap->n+1)*sizeof(MatScalar),&spdot);CHKERRQ(ierr); 770e2cac8caSJed Brown ierr = PetscMemzero(spdot,(A->cmap->n+1)*sizeof(MatScalar));CHKERRQ(ierr); 7716973769bSHong Zhang #endif 7726973769bSHong Zhang 77381d82fe4SHong Zhang /* clear old values in C */ 774aa1aec99SHong Zhang if (!c->a){ 775aa1aec99SHong Zhang ierr = PetscMalloc((ci[cm]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr); 776aa1aec99SHong Zhang c->a = ca; 777aa1aec99SHong Zhang c->free_a = PETSC_TRUE; 778aa1aec99SHong Zhang } else { 779aa1aec99SHong Zhang ca = c->a; 780aa1aec99SHong Zhang } 78181d82fe4SHong Zhang ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr); 7825df89d91SHong Zhang 7831fa1209cSHong Zhang for (i=0; i<cm; i++) { 78481d82fe4SHong Zhang anzi = ai[i+1] - ai[i]; 78581d82fe4SHong Zhang acol = aj + ai[i]; 7866973769bSHong Zhang aval = aa + ai[i]; 7871fa1209cSHong Zhang cnzi = ci[i+1] - ci[i]; 7881fa1209cSHong Zhang ccol = cj + ci[i]; 7896973769bSHong Zhang cval = ca + ci[i]; 7901fa1209cSHong Zhang for (j=0; j<cnzi; j++){ 79181d82fe4SHong Zhang brow = ccol[j]; 79281d82fe4SHong Zhang bnzj = bi[brow+1] - bi[brow]; 79381d82fe4SHong Zhang bcol = bj + bi[brow]; 7946973769bSHong Zhang bval = ba + bi[brow]; 7956973769bSHong Zhang 79681d82fe4SHong Zhang /* perform sparse inner-product c(i,j)=A[i,:]*B[j,:]^T */ 7976973769bSHong Zhang #if defined(USE_ARRAY) 7986973769bSHong Zhang /* put ba to spdot array */ 7996973769bSHong Zhang for (nextb=0; nextb<bnzj; nextb++) spdot[bcol[nextb]] = bval[nextb]; 8006973769bSHong Zhang /* c(i,j)=A[i,:]*B[j,:]^T */ 8016973769bSHong Zhang for (nexta=0; nexta<anzi; nexta++){ 8026973769bSHong Zhang cval[j] += spdot[acol[nexta]]*aval[nexta]; 8036973769bSHong Zhang } 8046973769bSHong Zhang /* zero spdot array */ 8056973769bSHong Zhang for (nextb=0; nextb<bnzj; nextb++) spdot[bcol[nextb]] = 0.0; 8066973769bSHong Zhang #else 80781d82fe4SHong Zhang nexta = 0; nextb = 0; 80881d82fe4SHong Zhang while (nexta<anzi && nextb<bnzj){ 80981d82fe4SHong Zhang while (acol[nexta] < bcol[nextb] && nexta < anzi) nexta++; 81081d82fe4SHong Zhang if (nexta == anzi) break; 81181d82fe4SHong Zhang while (acol[nexta] > bcol[nextb] && nextb < bnzj) nextb++; 81281d82fe4SHong Zhang if (nextb == bnzj) break; 81381d82fe4SHong Zhang if (acol[nexta] == bcol[nextb]){ 8146973769bSHong Zhang cval[j] += aval[nexta]*bval[nextb]; 81581d82fe4SHong Zhang nexta++; nextb++; 81681d82fe4SHong Zhang flops += 2; 8171fa1209cSHong Zhang } 8181fa1209cSHong Zhang } 8196973769bSHong Zhang #endif 82081d82fe4SHong Zhang } 82181d82fe4SHong Zhang } 82281d82fe4SHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 82381d82fe4SHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 82481d82fe4SHong Zhang ierr = PetscLogFlops(flops);CHKERRQ(ierr); 8256973769bSHong Zhang #if defined(USE_ARRAY) 8266973769bSHong Zhang ierr = PetscFree(spdot); 8276973769bSHong Zhang #endif 8285df89d91SHong Zhang PetscFunctionReturn(0); 8295df89d91SHong Zhang } 8305df89d91SHong Zhang 8315df89d91SHong Zhang #undef __FUNCT__ 83275648e8dSHong Zhang #define __FUNCT__ "MatTransposeMatMult_SeqAIJ_SeqAIJ" 83375648e8dSHong Zhang PetscErrorCode MatTransposeMatMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) { 8345df89d91SHong Zhang PetscErrorCode ierr; 8355df89d91SHong Zhang 8365df89d91SHong Zhang PetscFunctionBegin; 8375df89d91SHong Zhang if (scall == MAT_INITIAL_MATRIX){ 83875648e8dSHong Zhang ierr = MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr); 8395df89d91SHong Zhang } 84075648e8dSHong Zhang ierr = MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr); 8415df89d91SHong Zhang PetscFunctionReturn(0); 8425df89d91SHong Zhang } 8435df89d91SHong Zhang 8445df89d91SHong Zhang #undef __FUNCT__ 84575648e8dSHong Zhang #define __FUNCT__ "MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ" 84675648e8dSHong Zhang PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C) 8475df89d91SHong Zhang { 848bc011b1eSHong Zhang PetscErrorCode ierr; 849bc011b1eSHong Zhang Mat At; 85038baddfdSBarry Smith PetscInt *ati,*atj; 851bc011b1eSHong Zhang 852bc011b1eSHong Zhang PetscFunctionBegin; 853bc011b1eSHong Zhang /* create symbolic At */ 854bc011b1eSHong Zhang ierr = MatGetSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr); 855d0f46423SBarry Smith ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,A->cmap->n,A->rmap->n,ati,atj,PETSC_NULL,&At);CHKERRQ(ierr); 856*a2f3521dSMark F. Adams At->rmap->bs = A->cmap->bs; 857*a2f3521dSMark F. Adams At->cmap->bs = B->cmap->bs; 858bc011b1eSHong Zhang 859bc011b1eSHong Zhang /* get symbolic C=At*B */ 860bc011b1eSHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(At,B,fill,C);CHKERRQ(ierr); 861bc011b1eSHong Zhang 862bc011b1eSHong Zhang /* clean up */ 8636bf464f9SBarry Smith ierr = MatDestroy(&At);CHKERRQ(ierr); 864bc011b1eSHong Zhang ierr = MatRestoreSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr); 865bc011b1eSHong Zhang PetscFunctionReturn(0); 866bc011b1eSHong Zhang } 867bc011b1eSHong Zhang 868bc011b1eSHong Zhang #undef __FUNCT__ 86975648e8dSHong Zhang #define __FUNCT__ "MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ" 87075648e8dSHong Zhang PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C) 871bc011b1eSHong Zhang { 872bc011b1eSHong Zhang PetscErrorCode ierr; 8730fbc74f4SHong Zhang Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data; 874d0f46423SBarry Smith PetscInt am=A->rmap->n,anzi,*ai=a->i,*aj=a->j,*bi=b->i,*bj,bnzi,nextb; 875d0f46423SBarry Smith PetscInt cm=C->rmap->n,*ci=c->i,*cj=c->j,crow,*cjj,i,j,k; 876d13dce4bSSatish Balay PetscLogDouble flops=0.0; 877aa1aec99SHong Zhang MatScalar *aa=a->a,*ba,*ca,*caj; 878bc011b1eSHong Zhang 879bc011b1eSHong Zhang PetscFunctionBegin; 880aa1aec99SHong Zhang if (!c->a){ 881aa1aec99SHong Zhang ierr = PetscMalloc((ci[cm]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr); 882aa1aec99SHong Zhang c->a = ca; 883aa1aec99SHong Zhang c->free_a = PETSC_TRUE; 884aa1aec99SHong Zhang } else { 885aa1aec99SHong Zhang ca = c->a; 886aa1aec99SHong Zhang } 887bc011b1eSHong Zhang /* clear old values in C */ 888bc011b1eSHong Zhang ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr); 889bc011b1eSHong Zhang 890bc011b1eSHong Zhang /* compute A^T*B using outer product (A^T)[:,i]*B[i,:] */ 891bc011b1eSHong Zhang for (i=0;i<am;i++) { 892bc011b1eSHong Zhang bj = b->j + bi[i]; 893bc011b1eSHong Zhang ba = b->a + bi[i]; 894bc011b1eSHong Zhang bnzi = bi[i+1] - bi[i]; 895bc011b1eSHong Zhang anzi = ai[i+1] - ai[i]; 896bc011b1eSHong Zhang for (j=0; j<anzi; j++) { 897bc011b1eSHong Zhang nextb = 0; 8980fbc74f4SHong Zhang crow = *aj++; 899bc011b1eSHong Zhang cjj = cj + ci[crow]; 900bc011b1eSHong Zhang caj = ca + ci[crow]; 901bc011b1eSHong Zhang /* perform sparse axpy operation. Note cjj includes bj. */ 902bc011b1eSHong Zhang for (k=0; nextb<bnzi; k++) { 9030fbc74f4SHong Zhang if (cjj[k] == *(bj+nextb)) { /* ccol == bcol */ 9040fbc74f4SHong Zhang caj[k] += (*aa)*(*(ba+nextb)); 905bc011b1eSHong Zhang nextb++; 906bc011b1eSHong Zhang } 907bc011b1eSHong Zhang } 908bc011b1eSHong Zhang flops += 2*bnzi; 9090fbc74f4SHong Zhang aa++; 910bc011b1eSHong Zhang } 911bc011b1eSHong Zhang } 912bc011b1eSHong Zhang 913bc011b1eSHong Zhang /* Assemble the final matrix and clean up */ 914bc011b1eSHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 915bc011b1eSHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 916bc011b1eSHong Zhang ierr = PetscLogFlops(flops);CHKERRQ(ierr); 917bc011b1eSHong Zhang PetscFunctionReturn(0); 918bc011b1eSHong Zhang } 9199123193aSHong Zhang 920ec01deb9SMatthew Knepley EXTERN_C_BEGIN 9219123193aSHong Zhang #undef __FUNCT__ 9229123193aSHong Zhang #define __FUNCT__ "MatMatMult_SeqAIJ_SeqDense" 9239123193aSHong Zhang PetscErrorCode MatMatMult_SeqAIJ_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 9249123193aSHong Zhang { 9259123193aSHong Zhang PetscErrorCode ierr; 9269123193aSHong Zhang 9279123193aSHong Zhang PetscFunctionBegin; 9289123193aSHong Zhang if (scall == MAT_INITIAL_MATRIX){ 9299123193aSHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqDense(A,B,fill,C);CHKERRQ(ierr); 9309123193aSHong Zhang } 9319123193aSHong Zhang ierr = MatMatMultNumeric_SeqAIJ_SeqDense(A,B,*C);CHKERRQ(ierr); 9329123193aSHong Zhang PetscFunctionReturn(0); 9339123193aSHong Zhang } 934ec01deb9SMatthew Knepley EXTERN_C_END 935edd81eecSMatthew Knepley 9369123193aSHong Zhang #undef __FUNCT__ 9379123193aSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqDense" 9389123193aSHong Zhang PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C) 9399123193aSHong Zhang { 9409123193aSHong Zhang PetscErrorCode ierr; 9419123193aSHong Zhang 9429123193aSHong Zhang PetscFunctionBegin; 9435a586d82SBarry Smith ierr = MatMatMultSymbolic_SeqDense_SeqDense(A,B,0.0,C);CHKERRQ(ierr); 9448cdbd757SHong Zhang (*C)->ops->matmult = MatMatMult_SeqAIJ_SeqDense; 9459123193aSHong Zhang PetscFunctionReturn(0); 9469123193aSHong Zhang } 9479123193aSHong Zhang 9489123193aSHong Zhang #undef __FUNCT__ 9499123193aSHong Zhang #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqDense" 9509123193aSHong Zhang PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqDense(Mat A,Mat B,Mat C) 9519123193aSHong Zhang { 952f32d5d43SBarry Smith Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data; 9539123193aSHong Zhang PetscErrorCode ierr; 954dd6ea824SBarry Smith PetscScalar *b,*c,r1,r2,r3,r4,*b1,*b2,*b3,*b4; 955dd6ea824SBarry Smith MatScalar *aa; 956d0f46423SBarry Smith PetscInt cm=C->rmap->n, cn=B->cmap->n, bm=B->rmap->n, col, i,j,n,*aj, am = A->rmap->n; 957f32d5d43SBarry Smith PetscInt am2 = 2*am, am3 = 3*am, bm4 = 4*bm,colam; 9589123193aSHong Zhang 9599123193aSHong Zhang PetscFunctionBegin; 960f32d5d43SBarry Smith if (!cm || !cn) PetscFunctionReturn(0); 961e32f2f54SBarry 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); 962e32f2f54SBarry 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); 963e32f2f54SBarry 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); 964f32d5d43SBarry Smith ierr = MatGetArray(B,&b);CHKERRQ(ierr); 965f32d5d43SBarry Smith ierr = MatGetArray(C,&c);CHKERRQ(ierr); 966f32d5d43SBarry Smith b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm; 967f32d5d43SBarry Smith for (col=0; col<cn-4; col += 4){ /* over columns of C */ 968f32d5d43SBarry Smith colam = col*am; 969f32d5d43SBarry Smith for (i=0; i<am; i++) { /* over rows of C in those columns */ 970f32d5d43SBarry Smith r1 = r2 = r3 = r4 = 0.0; 971f32d5d43SBarry Smith n = a->i[i+1] - a->i[i]; 972f32d5d43SBarry Smith aj = a->j + a->i[i]; 973f32d5d43SBarry Smith aa = a->a + a->i[i]; 974f32d5d43SBarry Smith for (j=0; j<n; j++) { 975f32d5d43SBarry Smith r1 += (*aa)*b1[*aj]; 976f32d5d43SBarry Smith r2 += (*aa)*b2[*aj]; 977f32d5d43SBarry Smith r3 += (*aa)*b3[*aj]; 978f32d5d43SBarry Smith r4 += (*aa++)*b4[*aj++]; 9799123193aSHong Zhang } 980f32d5d43SBarry Smith c[colam + i] = r1; 981f32d5d43SBarry Smith c[colam + am + i] = r2; 982f32d5d43SBarry Smith c[colam + am2 + i] = r3; 983f32d5d43SBarry Smith c[colam + am3 + i] = r4; 984f32d5d43SBarry Smith } 985f32d5d43SBarry Smith b1 += bm4; 986f32d5d43SBarry Smith b2 += bm4; 987f32d5d43SBarry Smith b3 += bm4; 988f32d5d43SBarry Smith b4 += bm4; 989f32d5d43SBarry Smith } 990f32d5d43SBarry Smith for (;col<cn; col++){ /* over extra columns of C */ 991f32d5d43SBarry Smith for (i=0; i<am; i++) { /* over rows of C in those columns */ 992f32d5d43SBarry Smith r1 = 0.0; 993f32d5d43SBarry Smith n = a->i[i+1] - a->i[i]; 994f32d5d43SBarry Smith aj = a->j + a->i[i]; 995f32d5d43SBarry Smith aa = a->a + a->i[i]; 996f32d5d43SBarry Smith 997f32d5d43SBarry Smith for (j=0; j<n; j++) { 998f32d5d43SBarry Smith r1 += (*aa++)*b1[*aj++]; 999f32d5d43SBarry Smith } 1000f32d5d43SBarry Smith c[col*am + i] = r1; 1001f32d5d43SBarry Smith } 1002f32d5d43SBarry Smith b1 += bm; 1003f32d5d43SBarry Smith } 1004dc0b31edSSatish Balay ierr = PetscLogFlops(cn*(2.0*a->nz));CHKERRQ(ierr); 1005f32d5d43SBarry Smith ierr = MatRestoreArray(B,&b);CHKERRQ(ierr); 1006f32d5d43SBarry Smith ierr = MatRestoreArray(C,&c);CHKERRQ(ierr); 1007f32d5d43SBarry Smith ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1008f32d5d43SBarry Smith ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1009f32d5d43SBarry Smith PetscFunctionReturn(0); 1010f32d5d43SBarry Smith } 1011f32d5d43SBarry Smith 1012f32d5d43SBarry Smith /* 10134324174eSBarry Smith Note very similar to MatMult_SeqAIJ(), should generate both codes from same base 1014f32d5d43SBarry Smith */ 1015f32d5d43SBarry Smith #undef __FUNCT__ 1016f32d5d43SBarry Smith #define __FUNCT__ "MatMatMultNumericAdd_SeqAIJ_SeqDense" 1017f32d5d43SBarry Smith PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat A,Mat B,Mat C) 1018f32d5d43SBarry Smith { 1019f32d5d43SBarry Smith Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data; 1020f32d5d43SBarry Smith PetscErrorCode ierr; 1021dd6ea824SBarry Smith PetscScalar *b,*c,r1,r2,r3,r4,*b1,*b2,*b3,*b4; 1022dd6ea824SBarry Smith MatScalar *aa; 1023d0f46423SBarry 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; 10244324174eSBarry Smith PetscInt am2 = 2*am, am3 = 3*am, bm4 = 4*bm,colam,*ridx; 1025f32d5d43SBarry Smith 1026f32d5d43SBarry Smith PetscFunctionBegin; 1027f32d5d43SBarry Smith if (!cm || !cn) PetscFunctionReturn(0); 1028f32d5d43SBarry Smith ierr = MatGetArray(B,&b);CHKERRQ(ierr); 1029f32d5d43SBarry Smith ierr = MatGetArray(C,&c);CHKERRQ(ierr); 1030f32d5d43SBarry Smith b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm; 10314324174eSBarry Smith 10324324174eSBarry Smith if (a->compressedrow.use){ /* use compressed row format */ 10334324174eSBarry Smith for (col=0; col<cn-4; col += 4){ /* over columns of C */ 10344324174eSBarry Smith colam = col*am; 10354324174eSBarry Smith arm = a->compressedrow.nrows; 10364324174eSBarry Smith ii = a->compressedrow.i; 10374324174eSBarry Smith ridx = a->compressedrow.rindex; 10384324174eSBarry Smith for (i=0; i<arm; i++) { /* over rows of C in those columns */ 10394324174eSBarry Smith r1 = r2 = r3 = r4 = 0.0; 10404324174eSBarry Smith n = ii[i+1] - ii[i]; 10414324174eSBarry Smith aj = a->j + ii[i]; 10424324174eSBarry Smith aa = a->a + ii[i]; 10434324174eSBarry Smith for (j=0; j<n; j++) { 10444324174eSBarry Smith r1 += (*aa)*b1[*aj]; 10454324174eSBarry Smith r2 += (*aa)*b2[*aj]; 10464324174eSBarry Smith r3 += (*aa)*b3[*aj]; 10474324174eSBarry Smith r4 += (*aa++)*b4[*aj++]; 10484324174eSBarry Smith } 10494324174eSBarry Smith c[colam + ridx[i]] += r1; 10504324174eSBarry Smith c[colam + am + ridx[i]] += r2; 10514324174eSBarry Smith c[colam + am2 + ridx[i]] += r3; 10524324174eSBarry Smith c[colam + am3 + ridx[i]] += r4; 10534324174eSBarry Smith } 10544324174eSBarry Smith b1 += bm4; 10554324174eSBarry Smith b2 += bm4; 10564324174eSBarry Smith b3 += bm4; 10574324174eSBarry Smith b4 += bm4; 10584324174eSBarry Smith } 10594324174eSBarry Smith for (;col<cn; col++){ /* over extra columns of C */ 10604324174eSBarry Smith colam = col*am; 10614324174eSBarry Smith arm = a->compressedrow.nrows; 10624324174eSBarry Smith ii = a->compressedrow.i; 10634324174eSBarry Smith ridx = a->compressedrow.rindex; 10644324174eSBarry Smith for (i=0; i<arm; i++) { /* over rows of C in those columns */ 10654324174eSBarry Smith r1 = 0.0; 10664324174eSBarry Smith n = ii[i+1] - ii[i]; 10674324174eSBarry Smith aj = a->j + ii[i]; 10684324174eSBarry Smith aa = a->a + ii[i]; 10694324174eSBarry Smith 10704324174eSBarry Smith for (j=0; j<n; j++) { 10714324174eSBarry Smith r1 += (*aa++)*b1[*aj++]; 10724324174eSBarry Smith } 10734324174eSBarry Smith c[col*am + ridx[i]] += r1; 10744324174eSBarry Smith } 10754324174eSBarry Smith b1 += bm; 10764324174eSBarry Smith } 10774324174eSBarry Smith } else { 1078f32d5d43SBarry Smith for (col=0; col<cn-4; col += 4){ /* over columns of C */ 1079f32d5d43SBarry Smith colam = col*am; 1080f32d5d43SBarry Smith for (i=0; i<am; i++) { /* over rows of C in those columns */ 1081f32d5d43SBarry Smith r1 = r2 = r3 = r4 = 0.0; 1082f32d5d43SBarry Smith n = a->i[i+1] - a->i[i]; 1083f32d5d43SBarry Smith aj = a->j + a->i[i]; 1084f32d5d43SBarry Smith aa = a->a + a->i[i]; 1085f32d5d43SBarry Smith for (j=0; j<n; j++) { 1086f32d5d43SBarry Smith r1 += (*aa)*b1[*aj]; 1087f32d5d43SBarry Smith r2 += (*aa)*b2[*aj]; 1088f32d5d43SBarry Smith r3 += (*aa)*b3[*aj]; 1089f32d5d43SBarry Smith r4 += (*aa++)*b4[*aj++]; 1090f32d5d43SBarry Smith } 1091f32d5d43SBarry Smith c[colam + i] += r1; 1092f32d5d43SBarry Smith c[colam + am + i] += r2; 1093f32d5d43SBarry Smith c[colam + am2 + i] += r3; 1094f32d5d43SBarry Smith c[colam + am3 + i] += r4; 1095f32d5d43SBarry Smith } 1096f32d5d43SBarry Smith b1 += bm4; 1097f32d5d43SBarry Smith b2 += bm4; 1098f32d5d43SBarry Smith b3 += bm4; 1099f32d5d43SBarry Smith b4 += bm4; 1100f32d5d43SBarry Smith } 1101f32d5d43SBarry Smith for (;col<cn; col++){ /* over extra columns of C */ 1102f32d5d43SBarry Smith for (i=0; i<am; i++) { /* over rows of C in those columns */ 1103f32d5d43SBarry Smith r1 = 0.0; 1104f32d5d43SBarry Smith n = a->i[i+1] - a->i[i]; 1105f32d5d43SBarry Smith aj = a->j + a->i[i]; 1106f32d5d43SBarry Smith aa = a->a + a->i[i]; 1107f32d5d43SBarry Smith 1108f32d5d43SBarry Smith for (j=0; j<n; j++) { 1109f32d5d43SBarry Smith r1 += (*aa++)*b1[*aj++]; 1110f32d5d43SBarry Smith } 1111f32d5d43SBarry Smith c[col*am + i] += r1; 1112f32d5d43SBarry Smith } 1113f32d5d43SBarry Smith b1 += bm; 1114f32d5d43SBarry Smith } 11154324174eSBarry Smith } 1116dc0b31edSSatish Balay ierr = PetscLogFlops(cn*2.0*a->nz);CHKERRQ(ierr); 1117f32d5d43SBarry Smith ierr = MatRestoreArray(B,&b);CHKERRQ(ierr); 1118f32d5d43SBarry Smith ierr = MatRestoreArray(C,&c);CHKERRQ(ierr); 11199123193aSHong Zhang PetscFunctionReturn(0); 11209123193aSHong Zhang } 1121b1683b59SHong Zhang 1122b1683b59SHong Zhang #undef __FUNCT__ 1123b9af6bddSHong Zhang #define __FUNCT__ "MatTransColoringApplySpToDen_SeqAIJ" 1124b9af6bddSHong Zhang PetscErrorCode MatTransColoringApplySpToDen_SeqAIJ(MatTransposeColoring coloring,Mat B,Mat Btdense) 1125c8db22f6SHong Zhang { 1126c8db22f6SHong Zhang PetscErrorCode ierr; 11272f5f1f90SHong Zhang Mat_SeqAIJ *b = (Mat_SeqAIJ*)B->data; 11282f5f1f90SHong Zhang Mat_SeqDense *btdense = (Mat_SeqDense*)Btdense->data; 11292f5f1f90SHong Zhang PetscInt *bi=b->i,*bj=b->j; 11302f5f1f90SHong Zhang PetscInt m=Btdense->rmap->n,n=Btdense->cmap->n,j,k,l,col,anz,*btcol,brow,ncolumns; 11312f5f1f90SHong Zhang MatScalar *btval,*btval_den,*ba=b->a; 11322f5f1f90SHong Zhang PetscInt *columns=coloring->columns,*colorforcol=coloring->colorforcol,ncolors=coloring->ncolors; 1133c8db22f6SHong Zhang 1134c8db22f6SHong Zhang PetscFunctionBegin; 11352f5f1f90SHong Zhang btval_den=btdense->v; 11362f5f1f90SHong Zhang ierr = PetscMemzero(btval_den,(m*n)*sizeof(MatScalar));CHKERRQ(ierr); 11372f5f1f90SHong Zhang for (k=0; k<ncolors; k++) { 11382f5f1f90SHong Zhang ncolumns = coloring->ncolumns[k]; 11392f5f1f90SHong Zhang for (l=0; l<ncolumns; l++) { /* insert a row of B to a column of Btdense */ 1140d2853540SHong Zhang col = *(columns + colorforcol[k] + l); 11412f5f1f90SHong Zhang btcol = bj + bi[col]; 11422f5f1f90SHong Zhang btval = ba + bi[col]; 11432f5f1f90SHong Zhang anz = bi[col+1] - bi[col]; 1144d2853540SHong Zhang for (j=0; j<anz; j++){ 11452f5f1f90SHong Zhang brow = btcol[j]; 11462f5f1f90SHong Zhang btval_den[brow] = btval[j]; 1147c8db22f6SHong Zhang } 1148c8db22f6SHong Zhang } 11492f5f1f90SHong Zhang btval_den += m; 1150c8db22f6SHong Zhang } 1151c8db22f6SHong Zhang PetscFunctionReturn(0); 1152c8db22f6SHong Zhang } 1153c8db22f6SHong Zhang 1154c8db22f6SHong Zhang #undef __FUNCT__ 1155b9af6bddSHong Zhang #define __FUNCT__ "MatTransColoringApplyDenToSp_SeqAIJ" 1156b9af6bddSHong Zhang PetscErrorCode MatTransColoringApplyDenToSp_SeqAIJ(MatTransposeColoring matcoloring,Mat Cden,Mat Csp) 11578972f759SHong Zhang { 11588972f759SHong Zhang PetscErrorCode ierr; 1159b2d2b04fSHong Zhang Mat_SeqAIJ *csp = (Mat_SeqAIJ*)Csp->data; 11602f5f1f90SHong Zhang PetscInt k,l,*row,*idx,m,ncolors=matcoloring->ncolors,nrows; 1161b2d2b04fSHong Zhang PetscScalar *ca_den,*cp_den,*ca=csp->a; 11622f5f1f90SHong Zhang PetscInt *rows=matcoloring->rows,*spidx=matcoloring->columnsforspidx,*colorforrow=matcoloring->colorforrow; 11638972f759SHong Zhang 11648972f759SHong Zhang PetscFunctionBegin; 11658972f759SHong Zhang ierr = MatGetLocalSize(Csp,&m,PETSC_NULL);CHKERRQ(ierr); 1166b2d2b04fSHong Zhang ierr = MatGetArray(Cden,&ca_den);CHKERRQ(ierr); 1167b2d2b04fSHong Zhang cp_den = ca_den; 11682f5f1f90SHong Zhang for (k=0; k<ncolors; k++) { 11692f5f1f90SHong Zhang nrows = matcoloring->nrows[k]; 11702f5f1f90SHong Zhang row = rows + colorforrow[k]; 11712f5f1f90SHong Zhang idx = spidx + colorforrow[k]; 11722f5f1f90SHong Zhang for (l=0; l<nrows; l++){ 11732f5f1f90SHong Zhang ca[idx[l]] = cp_den[row[l]]; 1174b2d2b04fSHong Zhang } 1175b2d2b04fSHong Zhang cp_den += m; 1176b2d2b04fSHong Zhang } 1177b2d2b04fSHong Zhang ierr = MatRestoreArray(Cden,&ca_den);CHKERRQ(ierr); 11788972f759SHong Zhang PetscFunctionReturn(0); 11798972f759SHong Zhang } 11808972f759SHong Zhang 1181e99be685SHong Zhang /* 1182e99be685SHong Zhang MatGetColumnIJ_SeqAIJ_Color() and MatRestoreColumnIJ_SeqAIJ_Color() are customized from 1183e99be685SHong Zhang MatGetColumnIJ_SeqAIJ() and MatRestoreColumnIJ_SeqAIJ() by adding an output 1184e99be685SHong Zhang spidx[], index of a->j, to be used for setting 'columnsforspidx' in MatTransposeColoringCreate_SeqAIJ(). 1185e99be685SHong Zhang */ 1186e99be685SHong Zhang #undef __FUNCT__ 1187e99be685SHong Zhang #define __FUNCT__ "MatGetColumnIJ_SeqAIJ_Color" 1188e99be685SHong Zhang PetscErrorCode MatGetColumnIJ_SeqAIJ_Color(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool inodecompressed,PetscInt *nn,PetscInt *ia[],PetscInt *ja[],PetscInt *spidx[],PetscBool *done) 1189e99be685SHong Zhang { 1190e99be685SHong Zhang Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 1191e99be685SHong Zhang PetscErrorCode ierr; 1192e99be685SHong Zhang PetscInt i,*collengths,*cia,*cja,n = A->cmap->n,m = A->rmap->n; 1193e99be685SHong Zhang PetscInt nz = a->i[m],row,*jj,mr,col; 1194e99be685SHong Zhang PetscInt *cspidx; 1195e99be685SHong Zhang 1196e99be685SHong Zhang PetscFunctionBegin; 1197e99be685SHong Zhang *nn = n; 1198e99be685SHong Zhang if (!ia) PetscFunctionReturn(0); 1199e99be685SHong Zhang if (symmetric) { 1200e99be685SHong Zhang SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"MatGetColumnIJ_SeqAIJ_Color() not supported for the case symmetric"); 1201e99be685SHong Zhang ierr = MatToSymmetricIJ_SeqAIJ(A->rmap->n,a->i,a->j,0,oshift,ia,ja);CHKERRQ(ierr); 1202e99be685SHong Zhang } else { 1203e99be685SHong Zhang ierr = PetscMalloc((n+1)*sizeof(PetscInt),&collengths);CHKERRQ(ierr); 1204e99be685SHong Zhang ierr = PetscMemzero(collengths,n*sizeof(PetscInt));CHKERRQ(ierr); 1205e99be685SHong Zhang ierr = PetscMalloc((n+1)*sizeof(PetscInt),&cia);CHKERRQ(ierr); 1206e99be685SHong Zhang ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&cja);CHKERRQ(ierr); 1207e99be685SHong Zhang ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&cspidx);CHKERRQ(ierr); 1208e99be685SHong Zhang jj = a->j; 1209e99be685SHong Zhang for (i=0; i<nz; i++) { 1210e99be685SHong Zhang collengths[jj[i]]++; 1211e99be685SHong Zhang } 1212e99be685SHong Zhang cia[0] = oshift; 1213e99be685SHong Zhang for (i=0; i<n; i++) { 1214e99be685SHong Zhang cia[i+1] = cia[i] + collengths[i]; 1215e99be685SHong Zhang } 1216e99be685SHong Zhang ierr = PetscMemzero(collengths,n*sizeof(PetscInt));CHKERRQ(ierr); 1217e99be685SHong Zhang jj = a->j; 1218e99be685SHong Zhang for (row=0; row<m; row++) { 1219e99be685SHong Zhang mr = a->i[row+1] - a->i[row]; 1220e99be685SHong Zhang for (i=0; i<mr; i++) { 1221e99be685SHong Zhang col = *jj++; 1222e99be685SHong Zhang cspidx[cia[col] + collengths[col] - oshift] = a->i[row] + i; /* index of a->j */ 1223e99be685SHong Zhang cja[cia[col] + collengths[col]++ - oshift] = row + oshift; 1224e99be685SHong Zhang } 1225e99be685SHong Zhang } 1226e99be685SHong Zhang ierr = PetscFree(collengths);CHKERRQ(ierr); 1227e99be685SHong Zhang *ia = cia; *ja = cja; 1228e99be685SHong Zhang *spidx = cspidx; 1229e99be685SHong Zhang } 1230e99be685SHong Zhang PetscFunctionReturn(0); 1231e99be685SHong Zhang } 1232e99be685SHong Zhang 1233e99be685SHong Zhang #undef __FUNCT__ 1234e99be685SHong Zhang #define __FUNCT__ "MatRestoreColumnIJ_SeqAIJ_Color" 1235e99be685SHong Zhang PetscErrorCode MatRestoreColumnIJ_SeqAIJ_Color(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool inodecompressed,PetscInt *n,PetscInt *ia[],PetscInt *ja[],PetscInt *spidx[],PetscBool *done) 1236e99be685SHong Zhang { 1237e99be685SHong Zhang PetscErrorCode ierr; 1238e99be685SHong Zhang 1239e99be685SHong Zhang PetscFunctionBegin; 1240e99be685SHong Zhang if (!ia) PetscFunctionReturn(0); 1241e99be685SHong Zhang 1242e99be685SHong Zhang ierr = PetscFree(*ia);CHKERRQ(ierr); 1243e99be685SHong Zhang ierr = PetscFree(*ja);CHKERRQ(ierr); 1244e99be685SHong Zhang ierr = PetscFree(*spidx);CHKERRQ(ierr); 1245e99be685SHong Zhang PetscFunctionReturn(0); 1246e99be685SHong Zhang } 1247e99be685SHong Zhang 12488972f759SHong Zhang #undef __FUNCT__ 1249b9af6bddSHong Zhang #define __FUNCT__ "MatTransposeColoringCreate_SeqAIJ" 1250b9af6bddSHong Zhang PetscErrorCode MatTransposeColoringCreate_SeqAIJ(Mat mat,ISColoring iscoloring,MatTransposeColoring c) 1251b1683b59SHong Zhang { 1252b1683b59SHong Zhang PetscErrorCode ierr; 1253d2853540SHong Zhang PetscInt i,n,nrows,N,j,k,m,*row_idx,*ci,*cj,ncols,col,cm; 1254b1683b59SHong Zhang const PetscInt *is; 1255b28d80bdSHong Zhang PetscInt nis = iscoloring->n,*rowhit,bs = 1; 1256b1683b59SHong Zhang IS *isa; 1257b1683b59SHong Zhang PetscBool done; 1258b1683b59SHong Zhang PetscBool flg1,flg2; 1259b28d80bdSHong Zhang Mat_SeqAIJ *csp = (Mat_SeqAIJ*)mat->data; 1260e99be685SHong Zhang PetscInt *colorforrow,*rows,*rows_i,*columnsforspidx,*columnsforspidx_i,*idxhit,*spidx; 1261d2853540SHong Zhang PetscInt *colorforcol,*columns,*columns_i; 1262b1683b59SHong Zhang 1263b1683b59SHong Zhang PetscFunctionBegin; 1264b1683b59SHong Zhang ierr = ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);CHKERRQ(ierr); 1265e99be685SHong Zhang 1266b1683b59SHong Zhang /* this is ugly way to get blocksize but cannot call MatGetBlockSize() because AIJ can have bs > 1 */ 1267251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)mat,MATSEQBAIJ,&flg1);CHKERRQ(ierr); 1268251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)mat,MATMPIBAIJ,&flg2);CHKERRQ(ierr); 1269b1683b59SHong Zhang if (flg1 || flg2) { 1270b1683b59SHong Zhang ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr); 1271b1683b59SHong Zhang } 1272b1683b59SHong Zhang 1273b1683b59SHong Zhang N = mat->cmap->N/bs; 1274b1683b59SHong Zhang c->M = mat->rmap->N/bs; /* set total rows, columns and local rows */ 1275b1683b59SHong Zhang c->N = mat->cmap->N/bs; 1276b1683b59SHong Zhang c->m = mat->rmap->N/bs; 1277b1683b59SHong Zhang c->rstart = 0; 1278b1683b59SHong Zhang 1279b1683b59SHong Zhang c->ncolors = nis; 1280b1683b59SHong Zhang ierr = PetscMalloc(nis*sizeof(PetscInt),&c->ncolumns);CHKERRQ(ierr); 1281b28d80bdSHong Zhang ierr = PetscMalloc(nis*sizeof(PetscInt),&c->nrows);CHKERRQ(ierr); 1282d2853540SHong Zhang ierr = PetscMalloc2(csp->nz+1,PetscInt,&rows,csp->nz+1,PetscInt,&columnsforspidx);CHKERRQ(ierr); 1283d2853540SHong Zhang ierr = PetscMalloc((nis+1)*sizeof(PetscInt),&colorforrow);CHKERRQ(ierr); 1284d2853540SHong Zhang colorforrow[0] = 0; 1285d2853540SHong Zhang rows_i = rows; 1286d2853540SHong Zhang columnsforspidx_i = columnsforspidx; 1287b1683b59SHong Zhang 1288d2853540SHong Zhang ierr = PetscMalloc((nis+1)*sizeof(PetscInt),&colorforcol);CHKERRQ(ierr); 1289d2853540SHong Zhang ierr = PetscMalloc((N+1)*sizeof(PetscInt),&columns);CHKERRQ(ierr); 1290d2853540SHong Zhang colorforcol[0] = 0; 1291d2853540SHong Zhang columns_i = columns; 1292d2853540SHong Zhang 1293e99be685SHong Zhang ierr = MatGetColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,&done);CHKERRQ(ierr); /* column-wise storage of mat */ 1294b1683b59SHong Zhang if (!done) SETERRQ1(((PetscObject)mat)->comm,PETSC_ERR_SUP,"MatGetColumnIJ() not supported for matrix type %s",((PetscObject)mat)->type_name); 1295b1683b59SHong Zhang 1296b28d80bdSHong Zhang cm = c->m; 1297b28d80bdSHong Zhang ierr = PetscMalloc((cm+1)*sizeof(PetscInt),&rowhit);CHKERRQ(ierr); 1298e99be685SHong Zhang ierr = PetscMalloc((cm+1)*sizeof(PetscInt),&idxhit);CHKERRQ(ierr); 1299b1683b59SHong Zhang for (i=0; i<nis; i++) { 1300b1683b59SHong Zhang ierr = ISGetLocalSize(isa[i],&n);CHKERRQ(ierr); 1301b1683b59SHong Zhang ierr = ISGetIndices(isa[i],&is);CHKERRQ(ierr); 1302b1683b59SHong Zhang c->ncolumns[i] = n; 1303b1683b59SHong Zhang if (n) { 1304d2853540SHong Zhang ierr = PetscMemcpy(columns_i,is,n*sizeof(PetscInt));CHKERRQ(ierr); 1305b1683b59SHong Zhang } 1306d2853540SHong Zhang colorforcol[i+1] = colorforcol[i] + n; 1307d2853540SHong Zhang columns_i += n; 1308b1683b59SHong Zhang 1309b1683b59SHong Zhang /* fast, crude version requires O(N*N) work */ 1310e8354b3eSHong Zhang ierr = PetscMemzero(rowhit,cm*sizeof(PetscInt));CHKERRQ(ierr); 1311e99be685SHong Zhang 1312b1683b59SHong Zhang /* loop over columns*/ 1313b1683b59SHong Zhang for (j=0; j<n; j++) { 1314b1683b59SHong Zhang col = is[j]; 1315d2853540SHong Zhang row_idx = cj + ci[col]; 1316b1683b59SHong Zhang m = ci[col+1] - ci[col]; 1317b1683b59SHong Zhang /* loop over columns marking them in rowhit */ 1318b1683b59SHong Zhang for (k=0; k<m; k++) { 1319e99be685SHong Zhang idxhit[*row_idx] = spidx[ci[col] + k]; 1320d2853540SHong Zhang rowhit[*row_idx++] = col + 1; 1321b1683b59SHong Zhang } 1322b1683b59SHong Zhang } 1323b1683b59SHong Zhang /* count the number of hits */ 1324b1683b59SHong Zhang nrows = 0; 1325e8354b3eSHong Zhang for (j=0; j<cm; j++) { 1326b1683b59SHong Zhang if (rowhit[j]) nrows++; 1327b1683b59SHong Zhang } 1328b1683b59SHong Zhang c->nrows[i] = nrows; 1329d2853540SHong Zhang colorforrow[i+1] = colorforrow[i] + nrows; 1330d2853540SHong Zhang 1331b1683b59SHong Zhang nrows = 0; 1332e8354b3eSHong Zhang for (j=0; j<cm; j++) { 1333b1683b59SHong Zhang if (rowhit[j]) { 1334d2853540SHong Zhang rows_i[nrows] = j; 1335e99be685SHong Zhang columnsforspidx_i[nrows] = idxhit[j]; 1336b1683b59SHong Zhang nrows++; 1337b1683b59SHong Zhang } 1338b1683b59SHong Zhang } 1339b1683b59SHong Zhang ierr = ISRestoreIndices(isa[i],&is);CHKERRQ(ierr); 1340d2853540SHong Zhang rows_i += nrows; columnsforspidx_i += nrows; 1341b1683b59SHong Zhang } 1342e99be685SHong Zhang ierr = MatRestoreColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,&done);CHKERRQ(ierr); 1343b28d80bdSHong Zhang ierr = PetscFree(rowhit);CHKERRQ(ierr); 1344b1683b59SHong Zhang ierr = ISColoringRestoreIS(iscoloring,&isa);CHKERRQ(ierr); 1345d2853540SHong Zhang #if defined(PETSC_USE_DEBUG) 1346d2853540SHong Zhang if (csp->nz != colorforrow[nis]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"csp->nz %d != colorforrow[nis] %d",csp->nz,colorforrow[nis]); 1347d2853540SHong Zhang #endif 1348b28d80bdSHong Zhang 1349d2853540SHong Zhang c->colorforrow = colorforrow; 1350d2853540SHong Zhang c->rows = rows; 1351d2853540SHong Zhang c->columnsforspidx = columnsforspidx; 1352d2853540SHong Zhang c->colorforcol = colorforcol; 1353d2853540SHong Zhang c->columns = columns; 1354f94ccd6cSHong Zhang ierr = PetscFree(idxhit);CHKERRQ(ierr); 1355b1683b59SHong Zhang PetscFunctionReturn(0); 1356b1683b59SHong Zhang } 1357