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> 90ced3a2bSJed Brown #include <../src/mat/utils/petscheap.h> 10c6db04a5SJed Brown #include <petscbt.h> 11c6db04a5SJed Brown #include <../src/mat/impls/dense/seq/dense.h> /*I "petscmat.h" I*/ 127bab7c10SHong Zhang 136284ec50SHong Zhang #undef __FUNCT__ 145c66b693SKris Buschelman #define __FUNCT__ "MatMatMult_SeqAIJ_SeqAIJ" 1538baddfdSBarry Smith PetscErrorCode MatMatMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 1638baddfdSBarry Smith { 17dfbe8321SBarry Smith PetscErrorCode ierr; 188a07c6f1SJed Brown PetscBool scalable=PETSC_FALSE,scalable_fast=PETSC_FALSE,heap = PETSC_FALSE,btheap = PETSC_FALSE; 195c66b693SKris Buschelman 205c66b693SKris Buschelman PetscFunctionBegin; 2126be0446SHong Zhang if (scall == MAT_INITIAL_MATRIX){ 22152983bfSHong Zhang ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr); 23152983bfSHong Zhang ierr = PetscOptionsBool("-matmatmult_scalable","Use a scalable but slower C=A*B","",scalable,&scalable,PETSC_NULL);CHKERRQ(ierr); 243c50cad2SHong Zhang ierr = PetscOptionsBool("-matmatmult_scalable_fast","Use a scalable but slower C=A*B","",scalable_fast,&scalable_fast,PETSC_NULL);CHKERRQ(ierr); 250ced3a2bSJed Brown ierr = PetscOptionsBool("-matmatmult_heap","Use heap implementation of symbolic factorization C=A*B","",heap,&heap,PETSC_NULL);CHKERRQ(ierr); 268a07c6f1SJed Brown ierr = PetscOptionsBool("-matmatmult_btheap","Use btheap implementation of symbolic factorization C=A*B","",btheap,&btheap,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); 330ced3a2bSJed Brown } else if (heap) { 340ced3a2bSJed Brown ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_Heap(A,B,fill,C);CHKERRQ(ierr); 358a07c6f1SJed Brown } else if (btheap) { 368a07c6f1SJed Brown ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_BTHeap(A,B,fill,C);CHKERRQ(ierr); 3725023028SHong Zhang } else { 3826be0446SHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr); 3925023028SHong Zhang } 40598bc09dSHong Zhang ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr); 4126be0446SHong Zhang } 42598bc09dSHong Zhang ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr); 4301e47db4SHong Zhang ierr = (*(*C)->ops->matmultnumeric)(A,B,*C);CHKERRQ(ierr); 44598bc09dSHong Zhang ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr); 455c66b693SKris Buschelman PetscFunctionReturn(0); 465c66b693SKris Buschelman } 471c24bd37SHong Zhang 48e9e4536cSHong Zhang #undef __FUNCT__ 49b561aa9dSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ" 50b561aa9dSHong Zhang PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C) 51b561aa9dSHong Zhang { 52b561aa9dSHong Zhang PetscErrorCode ierr; 53b561aa9dSHong Zhang Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c; 54971236abSHong Zhang PetscInt *ai=a->i,*bi=b->i,*ci,*cj; 55c1ba5575SJed Brown PetscInt am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N; 56b561aa9dSHong Zhang PetscReal afill; 57fb908031SHong Zhang PetscInt i,j,anzi,brow,bnzj,cnzi,*bj,*aj,nlnk_max,*lnk,ndouble=0; 58fb908031SHong Zhang PetscBT lnkbt; 59fb908031SHong Zhang PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 60b561aa9dSHong Zhang 61b561aa9dSHong Zhang PetscFunctionBegin; 62bd958071SHong Zhang /* Get ci and cj */ 63bd958071SHong Zhang /*---------------*/ 64fb908031SHong Zhang /* Allocate ci array, arrays for fill computation and */ 65fb908031SHong Zhang /* free space for accumulating nonzero column info */ 66fb908031SHong Zhang ierr = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr); 67fb908031SHong Zhang ci[0] = 0; 68fb908031SHong Zhang 69fb908031SHong Zhang /* create and initialize a linked list */ 70fb908031SHong Zhang nlnk_max = a->rmax*b->rmax; 71fb908031SHong Zhang if (!nlnk_max || nlnk_max > bn) nlnk_max = bn; 72fb908031SHong Zhang ierr = PetscLLCondensedCreate(nlnk_max,bn,&lnk,&lnkbt);CHKERRQ(ierr); 73fb908031SHong Zhang 74fb908031SHong Zhang /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */ 75fb908031SHong Zhang ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr); 76fb908031SHong Zhang current_space = free_space; 77fb908031SHong Zhang 78fb908031SHong Zhang /* Determine ci and cj */ 79fb908031SHong Zhang for (i=0; i<am; i++) { 80fb908031SHong Zhang anzi = ai[i+1] - ai[i]; 81fb908031SHong Zhang aj = a->j + ai[i]; 82fb908031SHong Zhang for (j=0; j<anzi; j++){ 83fb908031SHong Zhang brow = aj[j]; 84fb908031SHong Zhang bnzj = bi[brow+1] - bi[brow]; 85fb908031SHong Zhang bj = b->j + bi[brow]; 86fb908031SHong Zhang /* add non-zero cols of B into the sorted linked list lnk */ 87fb908031SHong Zhang ierr = PetscLLCondensedAddSorted(bnzj,bj,lnk,lnkbt);CHKERRQ(ierr); 88fb908031SHong Zhang } 89fb908031SHong Zhang cnzi = lnk[0]; 90fb908031SHong Zhang 91fb908031SHong Zhang /* If free space is not available, make more free space */ 92fb908031SHong Zhang /* Double the amount of total space in the list */ 93fb908031SHong Zhang if (current_space->local_remaining<cnzi) { 94fb908031SHong Zhang ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 95fb908031SHong Zhang ndouble++; 96fb908031SHong Zhang } 97fb908031SHong Zhang 98fb908031SHong Zhang /* Copy data into free space, then initialize lnk */ 99fb908031SHong Zhang ierr = PetscLLCondensedClean(bn,cnzi,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 100fb908031SHong Zhang current_space->array += cnzi; 101fb908031SHong Zhang current_space->local_used += cnzi; 102fb908031SHong Zhang current_space->local_remaining -= cnzi; 103fb908031SHong Zhang ci[i+1] = ci[i] + cnzi; 104fb908031SHong Zhang } 105fb908031SHong Zhang 106fb908031SHong Zhang /* Column indices are in the list of free space */ 107fb908031SHong Zhang /* Allocate space for cj, initialize cj, and */ 108fb908031SHong Zhang /* destroy list of free space and other temporary array(s) */ 109fb908031SHong Zhang ierr = PetscMalloc((ci[am]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr); 110fb908031SHong Zhang ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr); 111fb908031SHong Zhang ierr = PetscLLCondensedDestroy(lnk,lnkbt);CHKERRQ(ierr); 112b561aa9dSHong Zhang 11326be0446SHong Zhang /* put together the new symbolic matrix */ 114aa1aec99SHong Zhang ierr = MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,am,bn,ci,cj,PETSC_NULL,C);CHKERRQ(ierr); 115a2f3521dSMark F. Adams (*C)->rmap->bs = A->rmap->bs; 116a2f3521dSMark F. Adams (*C)->cmap->bs = B->cmap->bs; 11758c24d83SHong Zhang 11858c24d83SHong Zhang /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 11958c24d83SHong Zhang /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */ 12058c24d83SHong Zhang c = (Mat_SeqAIJ *)((*C)->data); 121aa1aec99SHong Zhang c->free_a = PETSC_FALSE; 122e6b907acSBarry Smith c->free_ij = PETSC_TRUE; 12358c24d83SHong Zhang c->nonew = 0; 124002e8affSHong Zhang (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ; /* fast, needs non-scalable O(bn) array 'abdense' */ 1250b7e3e3dSHong Zhang 1267212da7cSHong Zhang /* set MatInfo */ 1277212da7cSHong Zhang afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5; 128f2b054eeSHong Zhang if (afill < 1.0) afill = 1.0; 1297212da7cSHong Zhang c->maxnz = ci[am]; 1307212da7cSHong Zhang c->nz = ci[am]; 131fb908031SHong Zhang (*C)->info.mallocs = ndouble; 1327212da7cSHong Zhang (*C)->info.fill_ratio_given = fill; 1337212da7cSHong Zhang (*C)->info.fill_ratio_needed = afill; 1347212da7cSHong Zhang 1357212da7cSHong Zhang #if defined(PETSC_USE_INFO) 1367212da7cSHong Zhang if (ci[am]) { 137fb908031SHong Zhang ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",ndouble,fill,afill);CHKERRQ(ierr); 138f2b054eeSHong Zhang ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr); 139f2b054eeSHong Zhang } else { 140f2b054eeSHong Zhang ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr); 141be0fcf8dSHong Zhang } 142f2b054eeSHong Zhang #endif 14358c24d83SHong Zhang PetscFunctionReturn(0); 14458c24d83SHong Zhang } 145d50806bdSBarry Smith 14626be0446SHong Zhang #undef __FUNCT__ 14726be0446SHong Zhang #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqAIJ" 148dfbe8321SBarry Smith PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C) 149d50806bdSBarry Smith { 150dfbe8321SBarry Smith PetscErrorCode ierr; 151d13dce4bSSatish Balay PetscLogDouble flops=0.0; 152d50806bdSBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 153d50806bdSBarry Smith Mat_SeqAIJ *b = (Mat_SeqAIJ *)B->data; 154d50806bdSBarry Smith Mat_SeqAIJ *c = (Mat_SeqAIJ *)C->data; 15538baddfdSBarry Smith PetscInt *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j; 156c58ca2e3SHong Zhang PetscInt am=A->rmap->n,cm=C->rmap->n; 157519eb980SHong Zhang PetscInt i,j,k,anzi,bnzi,cnzi,brow; 158aa1aec99SHong Zhang PetscScalar *aa=a->a,*ba=b->a,*baj,*ca,valtmp; 159aa1aec99SHong Zhang PetscScalar *ab_dense; 160d50806bdSBarry Smith 161d50806bdSBarry Smith PetscFunctionBegin; 162aa1aec99SHong Zhang /* printf("MatMatMultNumeric_SeqAIJ_SeqAIJ...ca %p\n",c->a); */ 163aa1aec99SHong Zhang if (!c->a){ /* first call of MatMatMultNumeric_SeqAIJ_SeqAIJ, allocate ca and matmult_abdense */ 164aa1aec99SHong Zhang ierr = PetscMalloc((ci[cm]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr); 165aa1aec99SHong Zhang c->a = ca; 166aa1aec99SHong Zhang c->free_a = PETSC_TRUE; 167aa1aec99SHong Zhang 168aa1aec99SHong Zhang ierr = PetscMalloc(B->cmap->N*sizeof(PetscScalar),&ab_dense);CHKERRQ(ierr); 169aa1aec99SHong Zhang ierr = PetscMemzero(ab_dense,B->cmap->N*sizeof(PetscScalar));CHKERRQ(ierr); 170aa1aec99SHong Zhang c->matmult_abdense = ab_dense; 171aa1aec99SHong Zhang } else { 172aa1aec99SHong Zhang ca = c->a; 173aa1aec99SHong Zhang ab_dense = c->matmult_abdense; 174aa1aec99SHong Zhang } 175aa1aec99SHong Zhang 176c124e916SHong Zhang /* clean old values in C */ 177c124e916SHong Zhang ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr); 178d50806bdSBarry Smith /* Traverse A row-wise. */ 179d50806bdSBarry Smith /* Build the ith row in C by summing over nonzero columns in A, */ 180d50806bdSBarry Smith /* the rows of B corresponding to nonzeros of A. */ 181d50806bdSBarry Smith for (i=0; i<am; i++) { 182d50806bdSBarry Smith anzi = ai[i+1] - ai[i]; 183d50806bdSBarry Smith for (j=0; j<anzi; j++) { 184519eb980SHong Zhang brow = aj[j]; 185d50806bdSBarry Smith bnzi = bi[brow+1] - bi[brow]; 186d50806bdSBarry Smith bjj = bj + bi[brow]; 187d50806bdSBarry Smith baj = ba + bi[brow]; 188519eb980SHong Zhang /* perform dense axpy */ 18936ec6d2dSHong Zhang valtmp = aa[j]; 190519eb980SHong Zhang for (k=0; k<bnzi; k++) { 19136ec6d2dSHong Zhang ab_dense[bjj[k]] += valtmp*baj[k]; 192519eb980SHong Zhang } 193519eb980SHong Zhang flops += 2*bnzi; 194519eb980SHong Zhang } 195c58ca2e3SHong Zhang aj += anzi; aa += anzi; 196c58ca2e3SHong Zhang 197c58ca2e3SHong Zhang cnzi = ci[i+1] - ci[i]; 198519eb980SHong Zhang for (k=0; k<cnzi; k++) { 199519eb980SHong Zhang ca[k] += ab_dense[cj[k]]; 200519eb980SHong Zhang ab_dense[cj[k]] = 0.0; /* zero ab_dense */ 201519eb980SHong Zhang } 202519eb980SHong Zhang flops += cnzi; 203519eb980SHong Zhang cj += cnzi; ca += cnzi; 204519eb980SHong Zhang } 205c58ca2e3SHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 206c58ca2e3SHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 207c58ca2e3SHong Zhang ierr = PetscLogFlops(flops);CHKERRQ(ierr); 208c58ca2e3SHong Zhang PetscFunctionReturn(0); 209c58ca2e3SHong Zhang } 210c58ca2e3SHong Zhang 211c58ca2e3SHong Zhang #undef __FUNCT__ 21225023028SHong Zhang #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable" 21325023028SHong Zhang PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(Mat A,Mat B,Mat C) 214c58ca2e3SHong Zhang { 215c58ca2e3SHong Zhang PetscErrorCode ierr; 216c58ca2e3SHong Zhang PetscLogDouble flops=0.0; 217c58ca2e3SHong Zhang Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 218c58ca2e3SHong Zhang Mat_SeqAIJ *b = (Mat_SeqAIJ *)B->data; 219c58ca2e3SHong Zhang Mat_SeqAIJ *c = (Mat_SeqAIJ *)C->data; 220c58ca2e3SHong Zhang PetscInt *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j; 221c58ca2e3SHong Zhang PetscInt am=A->rmap->N,cm=C->rmap->N; 222c58ca2e3SHong Zhang PetscInt i,j,k,anzi,bnzi,cnzi,brow; 22336ec6d2dSHong Zhang PetscScalar *aa=a->a,*ba=b->a,*baj,*ca=c->a,valtmp; 224c58ca2e3SHong Zhang PetscInt nextb; 225c58ca2e3SHong Zhang 226c58ca2e3SHong Zhang PetscFunctionBegin; 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); 330a2f3521dSMark F. Adams (*C)->rmap->bs = A->rmap->bs; 331a2f3521dSMark 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); 433a2f3521dSMark F. Adams (*C)->rmap->bs = A->rmap->bs; 434a2f3521dSMark 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 4640ced3a2bSJed Brown #undef __FUNCT__ 4650ced3a2bSJed Brown #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ_Heap" 4660ced3a2bSJed Brown PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Heap(Mat A,Mat B,PetscReal fill,Mat *C) 4670ced3a2bSJed Brown { 4680ced3a2bSJed Brown PetscErrorCode ierr; 4690ced3a2bSJed Brown Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c; 4700ced3a2bSJed Brown const PetscInt *ai=a->i,*bi=b->i,*aj=a->j,*bj=b->j; 4710ced3a2bSJed Brown PetscInt *ci,*cj,*bb; 4720ced3a2bSJed Brown PetscInt am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N; 4730ced3a2bSJed Brown PetscReal afill; 4740ced3a2bSJed Brown PetscInt i,j,col,ndouble = 0; 4750ced3a2bSJed Brown PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 4760ced3a2bSJed Brown PetscHeap h; 4770ced3a2bSJed Brown 4780ced3a2bSJed Brown PetscFunctionBegin; 4790ced3a2bSJed Brown /* Get ci and cj - same as MatMatMultSymbolic_SeqAIJ_SeqAIJ except using PetscLLxxx_Scalalbe() */ 4800ced3a2bSJed Brown /*---------------------------------------------------------------------------------------------*/ 4810ced3a2bSJed Brown /* Allocate arrays for fill computation and free space for accumulating nonzero column */ 4820ced3a2bSJed Brown ierr = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr); 4830ced3a2bSJed Brown ci[0] = 0; 4840ced3a2bSJed Brown 4850ced3a2bSJed Brown /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */ 4860ced3a2bSJed Brown ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr); 4870ced3a2bSJed Brown current_space = free_space; 4880ced3a2bSJed Brown 4890ced3a2bSJed Brown ierr = PetscHeapCreate(a->rmax,&h);CHKERRQ(ierr); 4900ced3a2bSJed Brown ierr = PetscMalloc(a->rmax*sizeof(PetscInt),&bb);CHKERRQ(ierr); 4910ced3a2bSJed Brown 4920ced3a2bSJed Brown /* Determine ci and cj */ 4930ced3a2bSJed Brown for (i=0; i<am; i++) { 4940ced3a2bSJed Brown const PetscInt anzi = ai[i+1] - ai[i]; /* number of nonzeros in this row of A, this is the number of rows of B that we merge */ 4950ced3a2bSJed Brown const PetscInt *acol = aj + ai[i]; /* column indices of nonzero entries in this row */ 4960ced3a2bSJed Brown ci[i+1] = ci[i]; 4970ced3a2bSJed Brown /* Populate the min heap */ 4980ced3a2bSJed Brown for (j=0; j<anzi; j++) { 4990ced3a2bSJed Brown bb[j] = bi[acol[j]]; /* bb points at the start of the row */ 5000ced3a2bSJed Brown if (bb[j] < bi[acol[j]+1]) { /* Add if row is nonempty */ 5010ced3a2bSJed Brown ierr = PetscHeapAdd(h,j,bj[bb[j]++]);CHKERRQ(ierr); 5020ced3a2bSJed Brown } 5030ced3a2bSJed Brown } 5040ced3a2bSJed Brown /* Pick off the min element, adding it to free space */ 5050ced3a2bSJed Brown ierr = PetscHeapPop(h,&j,&col);CHKERRQ(ierr); 5060ced3a2bSJed Brown while (j >= 0) { 5070ced3a2bSJed Brown if (current_space->local_remaining < 1) { /* double the size, but don't exceed 16 MiB */ 5080ced3a2bSJed Brown ierr = PetscFreeSpaceGet(PetscMin(2*current_space->total_array_size,16 << 20),¤t_space);CHKERRQ(ierr); 5090ced3a2bSJed Brown ndouble++; 5100ced3a2bSJed Brown } 5110ced3a2bSJed Brown *(current_space->array++) = col; 5120ced3a2bSJed Brown current_space->local_used++; 5130ced3a2bSJed Brown current_space->local_remaining--; 5140ced3a2bSJed Brown ci[i+1]++; 5150ced3a2bSJed Brown 5160ced3a2bSJed Brown /* stash if anything else remains in this row of B */ 5170ced3a2bSJed Brown if (bb[j] < bi[acol[j]+1]) {ierr = PetscHeapStash(h,j,bj[bb[j]++]);CHKERRQ(ierr);} 5180ced3a2bSJed Brown while (1) { /* pop and stash any other rows of B that also had an entry in this column */ 5190ced3a2bSJed Brown PetscInt j2,col2; 5200ced3a2bSJed Brown ierr = PetscHeapPeek(h,&j2,&col2);CHKERRQ(ierr); 5210ced3a2bSJed Brown if (col2 != col) break; 5220ced3a2bSJed Brown ierr = PetscHeapPop(h,&j2,&col2);CHKERRQ(ierr); 5230ced3a2bSJed Brown if (bb[j2] < bi[acol[j2]+1]) {ierr = PetscHeapStash(h,j2,bj[bb[j2]++]);CHKERRQ(ierr);} 5240ced3a2bSJed Brown } 5250ced3a2bSJed Brown /* Put any stashed elements back into the min heap */ 5260ced3a2bSJed Brown ierr = PetscHeapUnstash(h);CHKERRQ(ierr); 5270ced3a2bSJed Brown ierr = PetscHeapPop(h,&j,&col);CHKERRQ(ierr); 5280ced3a2bSJed Brown } 5290ced3a2bSJed Brown } 5300ced3a2bSJed Brown ierr = PetscFree(bb);CHKERRQ(ierr); 5310ced3a2bSJed Brown ierr = PetscHeapDestroy(&h);CHKERRQ(ierr); 5320ced3a2bSJed Brown 5330ced3a2bSJed Brown /* Column indices are in the list of free space */ 5340ced3a2bSJed Brown /* Allocate space for cj, initialize cj, and */ 5350ced3a2bSJed Brown /* destroy list of free space and other temporary array(s) */ 5360ced3a2bSJed Brown ierr = PetscMalloc(ci[am]*sizeof(PetscInt),&cj);CHKERRQ(ierr); 5370ced3a2bSJed Brown ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr); 5380ced3a2bSJed Brown 5390ced3a2bSJed Brown /* put together the new symbolic matrix */ 54089d95c1aSJed Brown ierr = MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,am,bn,ci,cj,PETSC_NULL,C);CHKERRQ(ierr); 5410ced3a2bSJed Brown (*C)->rmap->bs = A->rmap->bs; 5420ced3a2bSJed Brown (*C)->cmap->bs = B->cmap->bs; 5430ced3a2bSJed Brown 5440ced3a2bSJed Brown /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 5450ced3a2bSJed Brown /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */ 5460ced3a2bSJed Brown c = (Mat_SeqAIJ *)((*C)->data); 5470ced3a2bSJed Brown c->free_a = PETSC_TRUE; 5480ced3a2bSJed Brown c->free_ij = PETSC_TRUE; 5490ced3a2bSJed Brown c->nonew = 0; 55089d95c1aSJed Brown (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ; 5510ced3a2bSJed Brown 5520ced3a2bSJed Brown /* set MatInfo */ 5530ced3a2bSJed Brown afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5; 5540ced3a2bSJed Brown if (afill < 1.0) afill = 1.0; 5550ced3a2bSJed Brown c->maxnz = ci[am]; 5560ced3a2bSJed Brown c->nz = ci[am]; 5570ced3a2bSJed Brown (*C)->info.mallocs = ndouble; 5580ced3a2bSJed Brown (*C)->info.fill_ratio_given = fill; 5590ced3a2bSJed Brown (*C)->info.fill_ratio_needed = afill; 5600ced3a2bSJed Brown 5610ced3a2bSJed Brown #if defined(PETSC_USE_INFO) 5620ced3a2bSJed Brown if (ci[am]) { 5630ced3a2bSJed Brown ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",ndouble,fill,afill);CHKERRQ(ierr); 5640ced3a2bSJed Brown ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr); 5650ced3a2bSJed Brown } else { 5660ced3a2bSJed Brown ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr); 5670ced3a2bSJed Brown } 5680ced3a2bSJed Brown #endif 5690ced3a2bSJed Brown PetscFunctionReturn(0); 5700ced3a2bSJed Brown } 571e9e4536cSHong Zhang 5728a07c6f1SJed Brown #undef __FUNCT__ 5738a07c6f1SJed Brown #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ_BTHeap" 5748a07c6f1SJed Brown PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_BTHeap(Mat A,Mat B,PetscReal fill,Mat *C) 5758a07c6f1SJed Brown { 5768a07c6f1SJed Brown PetscErrorCode ierr; 5778a07c6f1SJed Brown Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c; 5788a07c6f1SJed Brown const PetscInt *ai=a->i,*bi=b->i,*aj=a->j,*bj=b->j; 5798a07c6f1SJed Brown PetscInt *ci,*cj,*bb; 5808a07c6f1SJed Brown PetscInt am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N; 5818a07c6f1SJed Brown PetscReal afill; 5828a07c6f1SJed Brown PetscInt i,j,col,ndouble = 0; 5838a07c6f1SJed Brown PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 5848a07c6f1SJed Brown PetscHeap h; 5858a07c6f1SJed Brown PetscBT bt; 5868a07c6f1SJed Brown 5878a07c6f1SJed Brown PetscFunctionBegin; 5888a07c6f1SJed Brown /* Get ci and cj - same as MatMatMultSymbolic_SeqAIJ_SeqAIJ except using PetscLLxxx_Scalalbe() */ 5898a07c6f1SJed Brown /*---------------------------------------------------------------------------------------------*/ 5908a07c6f1SJed Brown /* Allocate arrays for fill computation and free space for accumulating nonzero column */ 5918a07c6f1SJed Brown ierr = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr); 5928a07c6f1SJed Brown ci[0] = 0; 5938a07c6f1SJed Brown 5948a07c6f1SJed Brown /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */ 5958a07c6f1SJed Brown ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr); 5968a07c6f1SJed Brown current_space = free_space; 5978a07c6f1SJed Brown 5988a07c6f1SJed Brown ierr = PetscHeapCreate(a->rmax,&h);CHKERRQ(ierr); 5998a07c6f1SJed Brown ierr = PetscMalloc(a->rmax*sizeof(PetscInt),&bb);CHKERRQ(ierr); 6008a07c6f1SJed Brown ierr = PetscBTCreate(bn,&bt);CHKERRQ(ierr); 6018a07c6f1SJed Brown 6028a07c6f1SJed Brown /* Determine ci and cj */ 6038a07c6f1SJed Brown for (i=0; i<am; i++) { 6048a07c6f1SJed Brown const PetscInt anzi = ai[i+1] - ai[i]; /* number of nonzeros in this row of A, this is the number of rows of B that we merge */ 6058a07c6f1SJed Brown const PetscInt *acol = aj + ai[i]; /* column indices of nonzero entries in this row */ 6068a07c6f1SJed Brown const PetscInt *fptr = current_space->array; /* Save beginning of the row so we can clear the BT later */ 6078a07c6f1SJed Brown ci[i+1] = ci[i]; 6088a07c6f1SJed Brown /* Populate the min heap */ 6098a07c6f1SJed Brown for (j=0; j<anzi; j++) { 6108a07c6f1SJed Brown PetscInt brow = acol[j]; 6118a07c6f1SJed Brown for (bb[j] = bi[brow]; bb[j] < bi[brow+1]; bb[j]++) { 6128a07c6f1SJed Brown PetscInt bcol = bj[bb[j]]; 6138a07c6f1SJed Brown if (!PetscBTLookupSet(bt,bcol)) { /* new entry */ 6148a07c6f1SJed Brown ierr = PetscHeapAdd(h,j,bcol);CHKERRQ(ierr); 6158a07c6f1SJed Brown bb[j]++; 6168a07c6f1SJed Brown break; 6178a07c6f1SJed Brown } 6188a07c6f1SJed Brown } 6198a07c6f1SJed Brown } 6208a07c6f1SJed Brown /* Pick off the min element, adding it to free space */ 6218a07c6f1SJed Brown ierr = PetscHeapPop(h,&j,&col);CHKERRQ(ierr); 6228a07c6f1SJed Brown while (j >= 0) { 6238a07c6f1SJed Brown if (current_space->local_remaining < 1) { /* double the size, but don't exceed 16 MiB */ 6248a07c6f1SJed Brown fptr = PETSC_NULL; /* need PetscBTMemzero */ 6258a07c6f1SJed Brown ierr = PetscFreeSpaceGet(PetscMin(2*current_space->total_array_size,16 << 20),¤t_space);CHKERRQ(ierr); 6268a07c6f1SJed Brown ndouble++; 6278a07c6f1SJed Brown } 6288a07c6f1SJed Brown *(current_space->array++) = col; 6298a07c6f1SJed Brown current_space->local_used++; 6308a07c6f1SJed Brown current_space->local_remaining--; 6318a07c6f1SJed Brown ci[i+1]++; 6328a07c6f1SJed Brown 6338a07c6f1SJed Brown /* stash if anything else remains in this row of B */ 6348a07c6f1SJed Brown for ( ; bb[j] < bi[acol[j]+1]; bb[j]++) { 6358a07c6f1SJed Brown PetscInt bcol = bj[bb[j]]; 6368a07c6f1SJed Brown if (!PetscBTLookupSet(bt,bcol)) { /* new entry */ 6378a07c6f1SJed Brown ierr = PetscHeapAdd(h,j,bcol);CHKERRQ(ierr); 6388a07c6f1SJed Brown bb[j]++; 6398a07c6f1SJed Brown break; 6408a07c6f1SJed Brown } 6418a07c6f1SJed Brown } 6428a07c6f1SJed Brown ierr = PetscHeapPop(h,&j,&col);CHKERRQ(ierr); 6438a07c6f1SJed Brown } 6448a07c6f1SJed Brown if (fptr) { /* Clear the bits for this row */ 6458a07c6f1SJed Brown for ( ; fptr<current_space->array; fptr++) {ierr = PetscBTClear(bt,*fptr);CHKERRQ(ierr);} 6468a07c6f1SJed Brown } else { /* We reallocated so we don't remember (easily) how to clear only the bits we changed */ 6478a07c6f1SJed Brown ierr = PetscBTMemzero(bn,bt);CHKERRQ(ierr); 6488a07c6f1SJed Brown } 6498a07c6f1SJed Brown } 6508a07c6f1SJed Brown ierr = PetscFree(bb);CHKERRQ(ierr); 6518a07c6f1SJed Brown ierr = PetscHeapDestroy(&h);CHKERRQ(ierr); 6528a07c6f1SJed Brown ierr = PetscBTDestroy(&bt);CHKERRQ(ierr); 6538a07c6f1SJed Brown 6548a07c6f1SJed Brown /* Column indices are in the list of free space */ 6558a07c6f1SJed Brown /* Allocate space for cj, initialize cj, and */ 6568a07c6f1SJed Brown /* destroy list of free space and other temporary array(s) */ 6578a07c6f1SJed Brown ierr = PetscMalloc(ci[am]*sizeof(PetscInt),&cj);CHKERRQ(ierr); 6588a07c6f1SJed Brown ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr); 6598a07c6f1SJed Brown 6608a07c6f1SJed Brown /* put together the new symbolic matrix */ 66189d95c1aSJed Brown ierr = MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,am,bn,ci,cj,PETSC_NULL,C);CHKERRQ(ierr); 6628a07c6f1SJed Brown (*C)->rmap->bs = A->rmap->bs; 6638a07c6f1SJed Brown (*C)->cmap->bs = B->cmap->bs; 6648a07c6f1SJed Brown 6658a07c6f1SJed Brown /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 6668a07c6f1SJed Brown /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */ 6678a07c6f1SJed Brown c = (Mat_SeqAIJ *)((*C)->data); 6688a07c6f1SJed Brown c->free_a = PETSC_TRUE; 6698a07c6f1SJed Brown c->free_ij = PETSC_TRUE; 6708a07c6f1SJed Brown c->nonew = 0; 67189d95c1aSJed Brown (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ; 6728a07c6f1SJed Brown 6738a07c6f1SJed Brown /* set MatInfo */ 6748a07c6f1SJed Brown afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5; 6758a07c6f1SJed Brown if (afill < 1.0) afill = 1.0; 6768a07c6f1SJed Brown c->maxnz = ci[am]; 6778a07c6f1SJed Brown c->nz = ci[am]; 6788a07c6f1SJed Brown (*C)->info.mallocs = ndouble; 6798a07c6f1SJed Brown (*C)->info.fill_ratio_given = fill; 6808a07c6f1SJed Brown (*C)->info.fill_ratio_needed = afill; 6818a07c6f1SJed Brown 6828a07c6f1SJed Brown #if defined(PETSC_USE_INFO) 6838a07c6f1SJed Brown if (ci[am]) { 6848a07c6f1SJed Brown ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",ndouble,fill,afill);CHKERRQ(ierr); 6858a07c6f1SJed Brown ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr); 6868a07c6f1SJed Brown } else { 6878a07c6f1SJed Brown ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr); 6888a07c6f1SJed Brown } 6898a07c6f1SJed Brown #endif 6908a07c6f1SJed Brown PetscFunctionReturn(0); 6918a07c6f1SJed Brown } 6928a07c6f1SJed Brown 693d2853540SHong Zhang /* This routine is not used. Should be removed! */ 694bc011b1eSHong Zhang #undef __FUNCT__ 6956fc122caSHong Zhang #define __FUNCT__ "MatMatTransposeMult_SeqAIJ_SeqAIJ" 6966fc122caSHong Zhang PetscErrorCode MatMatTransposeMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 6975df89d91SHong Zhang { 698bc011b1eSHong Zhang PetscErrorCode ierr; 699bc011b1eSHong Zhang 700bc011b1eSHong Zhang PetscFunctionBegin; 701bc011b1eSHong Zhang if (scall == MAT_INITIAL_MATRIX){ 7026fc122caSHong Zhang ierr = MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr); 703bc011b1eSHong Zhang } 7046fc122caSHong Zhang ierr = MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr); 705bc011b1eSHong Zhang PetscFunctionReturn(0); 706bc011b1eSHong Zhang } 707bc011b1eSHong Zhang 708bc011b1eSHong Zhang #undef __FUNCT__ 709e286af10SHong Zhang #define __FUNCT__ "PetscContainerDestroy_Mat_MatMatTransMult" 710e286af10SHong Zhang PetscErrorCode PetscContainerDestroy_Mat_MatMatTransMult(void *ptr) 7112128a86cSHong Zhang { 7122128a86cSHong Zhang PetscErrorCode ierr; 713e286af10SHong Zhang Mat_MatMatTransMult *multtrans=(Mat_MatMatTransMult*)ptr; 7142128a86cSHong Zhang 7152128a86cSHong Zhang PetscFunctionBegin; 716b9af6bddSHong Zhang ierr = MatTransposeColoringDestroy(&multtrans->matcoloring);CHKERRQ(ierr); 7172128a86cSHong Zhang ierr = MatDestroy(&multtrans->Bt_den);CHKERRQ(ierr); 7182128a86cSHong Zhang ierr = MatDestroy(&multtrans->ABt_den);CHKERRQ(ierr); 7192128a86cSHong Zhang ierr = PetscFree(multtrans);CHKERRQ(ierr); 7202128a86cSHong Zhang PetscFunctionReturn(0); 7212128a86cSHong Zhang } 7222128a86cSHong Zhang 7232128a86cSHong Zhang #undef __FUNCT__ 7242128a86cSHong Zhang #define __FUNCT__ "MatDestroy_SeqAIJ_MatMatMultTrans" 7252128a86cSHong Zhang PetscErrorCode MatDestroy_SeqAIJ_MatMatMultTrans(Mat A) 7262128a86cSHong Zhang { 7272128a86cSHong Zhang PetscErrorCode ierr; 7282128a86cSHong Zhang PetscContainer container; 729e286af10SHong Zhang Mat_MatMatTransMult *multtrans=PETSC_NULL; 7302128a86cSHong Zhang 7312128a86cSHong Zhang PetscFunctionBegin; 732e286af10SHong Zhang ierr = PetscObjectQuery((PetscObject)A,"Mat_MatMatTransMult",(PetscObject *)&container);CHKERRQ(ierr); 7332128a86cSHong Zhang if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit"); 7342128a86cSHong Zhang ierr = PetscContainerGetPointer(container,(void **)&multtrans);CHKERRQ(ierr); 7352128a86cSHong Zhang A->ops->destroy = multtrans->destroy; 7362128a86cSHong Zhang if (A->ops->destroy) { 7372128a86cSHong Zhang ierr = (*A->ops->destroy)(A);CHKERRQ(ierr); 7382128a86cSHong Zhang } 739e286af10SHong Zhang ierr = PetscObjectCompose((PetscObject)A,"Mat_MatMatTransMult",0);CHKERRQ(ierr); 7402128a86cSHong Zhang PetscFunctionReturn(0); 7412128a86cSHong Zhang } 7422128a86cSHong Zhang 7432128a86cSHong Zhang #undef __FUNCT__ 7446fc122caSHong Zhang #define __FUNCT__ "MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ" 7456fc122caSHong Zhang PetscErrorCode MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C) 746bc011b1eSHong Zhang { 7475df89d91SHong Zhang PetscErrorCode ierr; 74881d82fe4SHong Zhang Mat Bt; 74981d82fe4SHong Zhang PetscInt *bti,*btj; 750e286af10SHong Zhang Mat_MatMatTransMult *multtrans; 7512128a86cSHong Zhang PetscContainer container; 752d2853540SHong Zhang 75381d82fe4SHong Zhang PetscFunctionBegin; 75481d82fe4SHong Zhang /* create symbolic Bt */ 75581d82fe4SHong Zhang ierr = MatGetSymbolicTranspose_SeqAIJ(B,&bti,&btj);CHKERRQ(ierr); 75681d82fe4SHong Zhang ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,B->cmap->n,B->rmap->n,bti,btj,PETSC_NULL,&Bt);CHKERRQ(ierr); 757c3e5699bSSatish Balay Bt->rmap->bs = A->cmap->bs; 758c3e5699bSSatish Balay Bt->cmap->bs = B->cmap->bs; 75981d82fe4SHong Zhang 76081d82fe4SHong Zhang /* get symbolic C=A*Bt */ 76181d82fe4SHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,Bt,fill,C);CHKERRQ(ierr); 76281d82fe4SHong Zhang 7632128a86cSHong Zhang /* create a supporting struct for reuse intermidiate dense matrices with matcoloring */ 764e286af10SHong Zhang ierr = PetscNew(Mat_MatMatTransMult,&multtrans);CHKERRQ(ierr); 7652128a86cSHong Zhang 7662128a86cSHong Zhang /* attach the supporting struct to C */ 7672128a86cSHong Zhang ierr = PetscContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr); 7682128a86cSHong Zhang ierr = PetscContainerSetPointer(container,multtrans);CHKERRQ(ierr); 769e286af10SHong Zhang ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_Mat_MatMatTransMult);CHKERRQ(ierr); 770e286af10SHong Zhang ierr = PetscObjectCompose((PetscObject)(*C),"Mat_MatMatTransMult",(PetscObject)container);CHKERRQ(ierr); 7712128a86cSHong Zhang ierr = PetscContainerDestroy(&container);CHKERRQ(ierr); 7722128a86cSHong Zhang 7732128a86cSHong Zhang multtrans->usecoloring = PETSC_FALSE; 7742128a86cSHong Zhang multtrans->destroy = (*C)->ops->destroy; 7752128a86cSHong Zhang (*C)->ops->destroy = MatDestroy_SeqAIJ_MatMatMultTrans; 7762128a86cSHong Zhang 777f400d4cbSHong Zhang ierr = PetscOptionsGetBool(PETSC_NULL,"-matmattransmult_color",&multtrans->usecoloring,PETSC_NULL);CHKERRQ(ierr); 7782128a86cSHong Zhang if (multtrans->usecoloring){ 779b9af6bddSHong Zhang /* Create MatTransposeColoring from symbolic C=A*B^T */ 780b9af6bddSHong Zhang MatTransposeColoring matcoloring; 7812128a86cSHong Zhang ISColoring iscoloring; 7822128a86cSHong Zhang Mat Bt_dense,C_dense; 783e8354b3eSHong Zhang 784502de53fSHong Zhang ierr = MatGetColoring(*C,MATCOLORINGLF,&iscoloring);CHKERRQ(ierr); 785b9af6bddSHong Zhang ierr = MatTransposeColoringCreate(*C,iscoloring,&matcoloring);CHKERRQ(ierr); 7862128a86cSHong Zhang multtrans->matcoloring = matcoloring; 7872128a86cSHong Zhang ierr = ISColoringDestroy(&iscoloring);CHKERRQ(ierr); 7882128a86cSHong Zhang 7892128a86cSHong Zhang /* Create Bt_dense and C_dense = A*Bt_dense */ 7902128a86cSHong Zhang ierr = MatCreate(PETSC_COMM_SELF,&Bt_dense);CHKERRQ(ierr); 7912128a86cSHong Zhang ierr = MatSetSizes(Bt_dense,A->cmap->n,matcoloring->ncolors,A->cmap->n,matcoloring->ncolors);CHKERRQ(ierr); 7922128a86cSHong Zhang ierr = MatSetType(Bt_dense,MATSEQDENSE);CHKERRQ(ierr); 7932128a86cSHong Zhang ierr = MatSeqDenseSetPreallocation(Bt_dense,PETSC_NULL);CHKERRQ(ierr); 7942128a86cSHong Zhang Bt_dense->assembled = PETSC_TRUE; 7952128a86cSHong Zhang multtrans->Bt_den = Bt_dense; 7962128a86cSHong Zhang 7972128a86cSHong Zhang ierr = MatCreate(PETSC_COMM_SELF,&C_dense);CHKERRQ(ierr); 7982128a86cSHong Zhang ierr = MatSetSizes(C_dense,A->rmap->n,matcoloring->ncolors,A->rmap->n,matcoloring->ncolors);CHKERRQ(ierr); 7992128a86cSHong Zhang ierr = MatSetType(C_dense,MATSEQDENSE);CHKERRQ(ierr); 8002128a86cSHong Zhang ierr = MatSeqDenseSetPreallocation(C_dense,PETSC_NULL);CHKERRQ(ierr); 8012128a86cSHong Zhang Bt_dense->assembled = PETSC_TRUE; 8022128a86cSHong Zhang multtrans->ABt_den = C_dense; 803f94ccd6cSHong Zhang 804f94ccd6cSHong Zhang #if defined(PETSC_USE_INFO) 8051ce71dffSSatish Balay { 806f94ccd6cSHong Zhang Mat_SeqAIJ *c=(Mat_SeqAIJ*)(*C)->data; 807f94ccd6cSHong 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)); 8081ce71dffSSatish Balay } 809f94ccd6cSHong Zhang #endif 8102128a86cSHong Zhang } 811e99be685SHong Zhang /* clean up */ 812e99be685SHong Zhang ierr = MatDestroy(&Bt);CHKERRQ(ierr); 813e99be685SHong Zhang ierr = MatRestoreSymbolicTranspose_SeqAIJ(B,&bti,&btj);CHKERRQ(ierr); 8142128a86cSHong Zhang 815f94ccd6cSHong Zhang 816f94ccd6cSHong Zhang 81781d82fe4SHong Zhang #if defined(INEFFICIENT_ALGORITHM) 81881d82fe4SHong Zhang /* The algorithm below computes am*bm sparse inner-product - inefficient! It will be deleted later. */ 8191fa1209cSHong Zhang PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 8201fa1209cSHong Zhang Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c; 8211fa1209cSHong Zhang PetscInt *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*ci,*cj,*acol,*bcol; 8221fa1209cSHong Zhang PetscInt am=A->rmap->N,bm=B->rmap->N; 8231fa1209cSHong Zhang PetscInt i,j,anzi,bnzj,cnzi,nlnk,*lnk,nspacedouble=0,ka,kb,index[1]; 8241fa1209cSHong Zhang MatScalar *ca; 8251fa1209cSHong Zhang PetscBT lnkbt; 8261fa1209cSHong Zhang PetscReal afill; 8275df89d91SHong Zhang 8281fa1209cSHong Zhang /* Allocate row pointer array ci */ 8291fa1209cSHong Zhang ierr = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr); 8301fa1209cSHong Zhang ci[0] = 0; 8311fa1209cSHong Zhang 8321fa1209cSHong Zhang /* Create and initialize a linked list for C columns */ 8331fa1209cSHong Zhang nlnk = bm+1; 8341fa1209cSHong Zhang ierr = PetscLLCreate(bm,bm,nlnk,lnk,lnkbt);CHKERRQ(ierr); 8351fa1209cSHong Zhang 8361fa1209cSHong Zhang /* Initial FreeSpace with size fill*(nnz(A)+nnz(B)) */ 8371fa1209cSHong Zhang ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr); 8381fa1209cSHong Zhang current_space = free_space; 8391fa1209cSHong Zhang 8401fa1209cSHong Zhang /* Determine symbolic info for each row of the product A*B^T: */ 8411fa1209cSHong Zhang for (i=0; i<am; i++) { 8421fa1209cSHong Zhang anzi = ai[i+1] - ai[i]; 8431fa1209cSHong Zhang cnzi = 0; 8441fa1209cSHong Zhang acol = aj + ai[i]; 8451fa1209cSHong Zhang for (j=0; j<bm; j++){ 8461fa1209cSHong Zhang bnzj = bi[j+1] - bi[j]; 8471fa1209cSHong Zhang bcol= bj + bi[j]; 84881d82fe4SHong Zhang /* sparse inner-product c(i,j)=A[i,:]*B[j,:]^T */ 8491fa1209cSHong Zhang ka = 0; kb = 0; 8501fa1209cSHong Zhang while (ka < anzi && kb < bnzj){ 85181d82fe4SHong Zhang while (acol[ka] < bcol[kb] && ka < anzi) ka++; 85281d82fe4SHong Zhang if (ka == anzi) break; 85381d82fe4SHong Zhang while (acol[ka] > bcol[kb] && kb < bnzj) kb++; 85481d82fe4SHong Zhang if (kb == bnzj) break; 8551fa1209cSHong Zhang if (acol[ka] == bcol[kb]){ /* add nonzero c(i,j) to lnk */ 8561fa1209cSHong Zhang index[0] = j; 8571fa1209cSHong Zhang ierr = PetscLLAdd(1,index,bm,nlnk,lnk,lnkbt);CHKERRQ(ierr); 8581fa1209cSHong Zhang cnzi++; 8591fa1209cSHong Zhang break; 8601fa1209cSHong Zhang } 8611fa1209cSHong Zhang } 8621fa1209cSHong Zhang } 8631fa1209cSHong Zhang 8641fa1209cSHong Zhang /* If free space is not available, make more free space */ 8651fa1209cSHong Zhang /* Double the amount of total space in the list */ 8661fa1209cSHong Zhang if (current_space->local_remaining<cnzi) { 8671fa1209cSHong Zhang ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 8681fa1209cSHong Zhang nspacedouble++; 8691fa1209cSHong Zhang } 8701fa1209cSHong Zhang 8711fa1209cSHong Zhang /* Copy data into free space, then initialize lnk */ 8721fa1209cSHong Zhang ierr = PetscLLClean(bm,bm,cnzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 8731fa1209cSHong Zhang current_space->array += cnzi; 8741fa1209cSHong Zhang current_space->local_used += cnzi; 8751fa1209cSHong Zhang current_space->local_remaining -= cnzi; 8761fa1209cSHong Zhang 8771fa1209cSHong Zhang ci[i+1] = ci[i] + cnzi; 8781fa1209cSHong Zhang } 8791fa1209cSHong Zhang 8801fa1209cSHong Zhang 8811fa1209cSHong Zhang /* Column indices are in the list of free space. 8821fa1209cSHong Zhang Allocate array cj, copy column indices to cj, and destroy list of free space */ 8831fa1209cSHong Zhang ierr = PetscMalloc((ci[am]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr); 8841fa1209cSHong Zhang ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr); 8851fa1209cSHong Zhang ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 8861fa1209cSHong Zhang 8871fa1209cSHong Zhang /* Allocate space for ca */ 8881fa1209cSHong Zhang ierr = PetscMalloc((ci[am]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr); 8891fa1209cSHong Zhang ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr); 8901fa1209cSHong Zhang 8911fa1209cSHong Zhang /* put together the new symbolic matrix */ 8921fa1209cSHong Zhang ierr = MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,am,bm,ci,cj,ca,C);CHKERRQ(ierr); 893a2f3521dSMark F. Adams (*C)->rmap->bs = A->cmap->bs; 894a2f3521dSMark F. Adams (*C)->cmap->bs = B->cmap->bs; 8951fa1209cSHong Zhang 8961fa1209cSHong Zhang /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 8971fa1209cSHong Zhang /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */ 8981fa1209cSHong Zhang c = (Mat_SeqAIJ *)((*C)->data); 8991fa1209cSHong Zhang c->free_a = PETSC_TRUE; 9001fa1209cSHong Zhang c->free_ij = PETSC_TRUE; 9011fa1209cSHong Zhang c->nonew = 0; 9021fa1209cSHong Zhang 9031fa1209cSHong Zhang /* set MatInfo */ 9041fa1209cSHong Zhang afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5; 9051fa1209cSHong Zhang if (afill < 1.0) afill = 1.0; 9061fa1209cSHong Zhang c->maxnz = ci[am]; 9071fa1209cSHong Zhang c->nz = ci[am]; 9081fa1209cSHong Zhang (*C)->info.mallocs = nspacedouble; 9091fa1209cSHong Zhang (*C)->info.fill_ratio_given = fill; 9101fa1209cSHong Zhang (*C)->info.fill_ratio_needed = afill; 9111fa1209cSHong Zhang 9121fa1209cSHong Zhang #if defined(PETSC_USE_INFO) 9131fa1209cSHong Zhang if (ci[am]) { 9141fa1209cSHong Zhang ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr); 9156fc122caSHong Zhang ierr = PetscInfo1((*C),"Use MatMatTransposeMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr); 9161fa1209cSHong Zhang } else { 9171fa1209cSHong Zhang ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr); 9181fa1209cSHong Zhang } 9191fa1209cSHong Zhang #endif 92081d82fe4SHong Zhang #endif 9215df89d91SHong Zhang PetscFunctionReturn(0); 9225df89d91SHong Zhang } 9235df89d91SHong Zhang 9246973769bSHong Zhang /* #define USE_ARRAY - for sparse dot product. Slower than !USE_ARRAY */ 9255df89d91SHong Zhang #undef __FUNCT__ 9266fc122caSHong Zhang #define __FUNCT__ "MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ" 9276fc122caSHong Zhang PetscErrorCode MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C) 9285df89d91SHong Zhang { 9295df89d91SHong Zhang PetscErrorCode ierr; 9305df89d91SHong Zhang Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data; 931e2cac8caSJed Brown PetscInt *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,anzi,bnzj,nexta,nextb,*acol,*bcol,brow; 93281d82fe4SHong Zhang PetscInt cm=C->rmap->n,*ci=c->i,*cj=c->j,i,j,cnzi,*ccol; 9335df89d91SHong Zhang PetscLogDouble flops=0.0; 934aa1aec99SHong Zhang MatScalar *aa=a->a,*aval,*ba=b->a,*bval,*ca,*cval; 935e286af10SHong Zhang Mat_MatMatTransMult *multtrans; 9362128a86cSHong Zhang PetscContainer container; 9376973769bSHong Zhang #if defined(USE_ARRAY) 9386973769bSHong Zhang MatScalar *spdot; 9396973769bSHong Zhang #endif 9405df89d91SHong Zhang 9415df89d91SHong Zhang PetscFunctionBegin; 942e286af10SHong Zhang ierr = PetscObjectQuery((PetscObject)C,"Mat_MatMatTransMult",(PetscObject *)&container);CHKERRQ(ierr); 9432128a86cSHong Zhang if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit"); 9442128a86cSHong Zhang ierr = PetscContainerGetPointer(container,(void **)&multtrans);CHKERRQ(ierr); 9452128a86cSHong Zhang if (multtrans->usecoloring){ 946b9af6bddSHong Zhang MatTransposeColoring matcoloring = multtrans->matcoloring; 947c8db22f6SHong Zhang Mat Bt_dense; 948c8db22f6SHong Zhang PetscInt m,n; 949a0b95ffbSSatish Balay Mat C_dense = multtrans->ABt_den; 950a0b95ffbSSatish Balay 9512128a86cSHong Zhang Bt_dense = multtrans->Bt_den; 952c8db22f6SHong Zhang ierr = MatGetLocalSize(Bt_dense,&m,&n);CHKERRQ(ierr); 953c8db22f6SHong Zhang 954b9af6bddSHong Zhang /* Get Bt_dense by Apply MatTransposeColoring to B */ 955b9af6bddSHong Zhang ierr = MatTransColoringApplySpToDen(matcoloring,B,Bt_dense);CHKERRQ(ierr); 956c8db22f6SHong Zhang 957c8db22f6SHong Zhang /* C_dense = A*Bt_dense */ 9582128a86cSHong Zhang ierr = MatMatMultNumeric_SeqAIJ_SeqDense(A,Bt_dense,C_dense);CHKERRQ(ierr); 959c8db22f6SHong Zhang 9602128a86cSHong Zhang /* Recover C from C_dense */ 961b9af6bddSHong Zhang ierr = MatTransColoringApplyDenToSp(matcoloring,C_dense,C);CHKERRQ(ierr); 962c8db22f6SHong Zhang PetscFunctionReturn(0); 963c8db22f6SHong Zhang } 964c8db22f6SHong Zhang 9656973769bSHong Zhang #if defined(USE_ARRAY) 9666973769bSHong Zhang /* allocate an array for implementing sparse inner-product */ 967e2cac8caSJed Brown ierr = PetscMalloc((A->cmap->n+1)*sizeof(MatScalar),&spdot);CHKERRQ(ierr); 968e2cac8caSJed Brown ierr = PetscMemzero(spdot,(A->cmap->n+1)*sizeof(MatScalar));CHKERRQ(ierr); 9696973769bSHong Zhang #endif 9706973769bSHong Zhang 97181d82fe4SHong Zhang /* clear old values in C */ 972aa1aec99SHong Zhang if (!c->a){ 973aa1aec99SHong Zhang ierr = PetscMalloc((ci[cm]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr); 974aa1aec99SHong Zhang c->a = ca; 975aa1aec99SHong Zhang c->free_a = PETSC_TRUE; 976aa1aec99SHong Zhang } else { 977aa1aec99SHong Zhang ca = c->a; 978aa1aec99SHong Zhang } 97981d82fe4SHong Zhang ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr); 9805df89d91SHong Zhang 9811fa1209cSHong Zhang for (i=0; i<cm; i++) { 98281d82fe4SHong Zhang anzi = ai[i+1] - ai[i]; 98381d82fe4SHong Zhang acol = aj + ai[i]; 9846973769bSHong Zhang aval = aa + ai[i]; 9851fa1209cSHong Zhang cnzi = ci[i+1] - ci[i]; 9861fa1209cSHong Zhang ccol = cj + ci[i]; 9876973769bSHong Zhang cval = ca + ci[i]; 9881fa1209cSHong Zhang for (j=0; j<cnzi; j++){ 98981d82fe4SHong Zhang brow = ccol[j]; 99081d82fe4SHong Zhang bnzj = bi[brow+1] - bi[brow]; 99181d82fe4SHong Zhang bcol = bj + bi[brow]; 9926973769bSHong Zhang bval = ba + bi[brow]; 9936973769bSHong Zhang 99481d82fe4SHong Zhang /* perform sparse inner-product c(i,j)=A[i,:]*B[j,:]^T */ 9956973769bSHong Zhang #if defined(USE_ARRAY) 9966973769bSHong Zhang /* put ba to spdot array */ 9976973769bSHong Zhang for (nextb=0; nextb<bnzj; nextb++) spdot[bcol[nextb]] = bval[nextb]; 9986973769bSHong Zhang /* c(i,j)=A[i,:]*B[j,:]^T */ 9996973769bSHong Zhang for (nexta=0; nexta<anzi; nexta++){ 10006973769bSHong Zhang cval[j] += spdot[acol[nexta]]*aval[nexta]; 10016973769bSHong Zhang } 10026973769bSHong Zhang /* zero spdot array */ 10036973769bSHong Zhang for (nextb=0; nextb<bnzj; nextb++) spdot[bcol[nextb]] = 0.0; 10046973769bSHong Zhang #else 100581d82fe4SHong Zhang nexta = 0; nextb = 0; 100681d82fe4SHong Zhang while (nexta<anzi && nextb<bnzj){ 100781d82fe4SHong Zhang while (acol[nexta] < bcol[nextb] && nexta < anzi) nexta++; 100881d82fe4SHong Zhang if (nexta == anzi) break; 100981d82fe4SHong Zhang while (acol[nexta] > bcol[nextb] && nextb < bnzj) nextb++; 101081d82fe4SHong Zhang if (nextb == bnzj) break; 101181d82fe4SHong Zhang if (acol[nexta] == bcol[nextb]){ 10126973769bSHong Zhang cval[j] += aval[nexta]*bval[nextb]; 101381d82fe4SHong Zhang nexta++; nextb++; 101481d82fe4SHong Zhang flops += 2; 10151fa1209cSHong Zhang } 10161fa1209cSHong Zhang } 10176973769bSHong Zhang #endif 101881d82fe4SHong Zhang } 101981d82fe4SHong Zhang } 102081d82fe4SHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 102181d82fe4SHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 102281d82fe4SHong Zhang ierr = PetscLogFlops(flops);CHKERRQ(ierr); 10236973769bSHong Zhang #if defined(USE_ARRAY) 10246973769bSHong Zhang ierr = PetscFree(spdot); 10256973769bSHong Zhang #endif 10265df89d91SHong Zhang PetscFunctionReturn(0); 10275df89d91SHong Zhang } 10285df89d91SHong Zhang 10295df89d91SHong Zhang #undef __FUNCT__ 103075648e8dSHong Zhang #define __FUNCT__ "MatTransposeMatMult_SeqAIJ_SeqAIJ" 103175648e8dSHong Zhang PetscErrorCode MatTransposeMatMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) { 10325df89d91SHong Zhang PetscErrorCode ierr; 10335df89d91SHong Zhang 10345df89d91SHong Zhang PetscFunctionBegin; 10355df89d91SHong Zhang if (scall == MAT_INITIAL_MATRIX){ 103675648e8dSHong Zhang ierr = MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr); 10375df89d91SHong Zhang } 103875648e8dSHong Zhang ierr = MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr); 10395df89d91SHong Zhang PetscFunctionReturn(0); 10405df89d91SHong Zhang } 10415df89d91SHong Zhang 10425df89d91SHong Zhang #undef __FUNCT__ 104375648e8dSHong Zhang #define __FUNCT__ "MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ" 104475648e8dSHong Zhang PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C) 10455df89d91SHong Zhang { 1046bc011b1eSHong Zhang PetscErrorCode ierr; 1047bc011b1eSHong Zhang Mat At; 104838baddfdSBarry Smith PetscInt *ati,*atj; 1049bc011b1eSHong Zhang 1050bc011b1eSHong Zhang PetscFunctionBegin; 1051bc011b1eSHong Zhang /* create symbolic At */ 1052bc011b1eSHong Zhang ierr = MatGetSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr); 1053d0f46423SBarry Smith ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,A->cmap->n,A->rmap->n,ati,atj,PETSC_NULL,&At);CHKERRQ(ierr); 1054a2f3521dSMark F. Adams At->rmap->bs = A->cmap->bs; 1055a2f3521dSMark F. Adams At->cmap->bs = B->cmap->bs; 1056bc011b1eSHong Zhang 1057bc011b1eSHong Zhang /* get symbolic C=At*B */ 1058bc011b1eSHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(At,B,fill,C);CHKERRQ(ierr); 1059bc011b1eSHong Zhang 1060bc011b1eSHong Zhang /* clean up */ 10616bf464f9SBarry Smith ierr = MatDestroy(&At);CHKERRQ(ierr); 1062bc011b1eSHong Zhang ierr = MatRestoreSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr); 1063bc011b1eSHong Zhang PetscFunctionReturn(0); 1064bc011b1eSHong Zhang } 1065bc011b1eSHong Zhang 1066bc011b1eSHong Zhang #undef __FUNCT__ 106775648e8dSHong Zhang #define __FUNCT__ "MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ" 106875648e8dSHong Zhang PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C) 1069bc011b1eSHong Zhang { 1070bc011b1eSHong Zhang PetscErrorCode ierr; 10710fbc74f4SHong Zhang Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data; 1072d0f46423SBarry Smith PetscInt am=A->rmap->n,anzi,*ai=a->i,*aj=a->j,*bi=b->i,*bj,bnzi,nextb; 1073d0f46423SBarry Smith PetscInt cm=C->rmap->n,*ci=c->i,*cj=c->j,crow,*cjj,i,j,k; 1074d13dce4bSSatish Balay PetscLogDouble flops=0.0; 1075aa1aec99SHong Zhang MatScalar *aa=a->a,*ba,*ca,*caj; 1076bc011b1eSHong Zhang 1077bc011b1eSHong Zhang PetscFunctionBegin; 1078aa1aec99SHong Zhang if (!c->a){ 1079aa1aec99SHong Zhang ierr = PetscMalloc((ci[cm]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr); 1080aa1aec99SHong Zhang c->a = ca; 1081aa1aec99SHong Zhang c->free_a = PETSC_TRUE; 1082aa1aec99SHong Zhang } else { 1083aa1aec99SHong Zhang ca = c->a; 1084aa1aec99SHong Zhang } 1085bc011b1eSHong Zhang /* clear old values in C */ 1086bc011b1eSHong Zhang ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr); 1087bc011b1eSHong Zhang 1088bc011b1eSHong Zhang /* compute A^T*B using outer product (A^T)[:,i]*B[i,:] */ 1089bc011b1eSHong Zhang for (i=0;i<am;i++) { 1090bc011b1eSHong Zhang bj = b->j + bi[i]; 1091bc011b1eSHong Zhang ba = b->a + bi[i]; 1092bc011b1eSHong Zhang bnzi = bi[i+1] - bi[i]; 1093bc011b1eSHong Zhang anzi = ai[i+1] - ai[i]; 1094bc011b1eSHong Zhang for (j=0; j<anzi; j++) { 1095bc011b1eSHong Zhang nextb = 0; 10960fbc74f4SHong Zhang crow = *aj++; 1097bc011b1eSHong Zhang cjj = cj + ci[crow]; 1098bc011b1eSHong Zhang caj = ca + ci[crow]; 1099bc011b1eSHong Zhang /* perform sparse axpy operation. Note cjj includes bj. */ 1100bc011b1eSHong Zhang for (k=0; nextb<bnzi; k++) { 11010fbc74f4SHong Zhang if (cjj[k] == *(bj+nextb)) { /* ccol == bcol */ 11020fbc74f4SHong Zhang caj[k] += (*aa)*(*(ba+nextb)); 1103bc011b1eSHong Zhang nextb++; 1104bc011b1eSHong Zhang } 1105bc011b1eSHong Zhang } 1106bc011b1eSHong Zhang flops += 2*bnzi; 11070fbc74f4SHong Zhang aa++; 1108bc011b1eSHong Zhang } 1109bc011b1eSHong Zhang } 1110bc011b1eSHong Zhang 1111bc011b1eSHong Zhang /* Assemble the final matrix and clean up */ 1112bc011b1eSHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1113bc011b1eSHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1114bc011b1eSHong Zhang ierr = PetscLogFlops(flops);CHKERRQ(ierr); 1115bc011b1eSHong Zhang PetscFunctionReturn(0); 1116bc011b1eSHong Zhang } 11179123193aSHong Zhang 1118ec01deb9SMatthew Knepley EXTERN_C_BEGIN 11199123193aSHong Zhang #undef __FUNCT__ 11209123193aSHong Zhang #define __FUNCT__ "MatMatMult_SeqAIJ_SeqDense" 11219123193aSHong Zhang PetscErrorCode MatMatMult_SeqAIJ_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 11229123193aSHong Zhang { 11239123193aSHong Zhang PetscErrorCode ierr; 11249123193aSHong Zhang 11259123193aSHong Zhang PetscFunctionBegin; 11269123193aSHong Zhang if (scall == MAT_INITIAL_MATRIX){ 11279123193aSHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqDense(A,B,fill,C);CHKERRQ(ierr); 11289123193aSHong Zhang } 11299123193aSHong Zhang ierr = MatMatMultNumeric_SeqAIJ_SeqDense(A,B,*C);CHKERRQ(ierr); 11309123193aSHong Zhang PetscFunctionReturn(0); 11319123193aSHong Zhang } 1132ec01deb9SMatthew Knepley EXTERN_C_END 1133edd81eecSMatthew Knepley 11349123193aSHong Zhang #undef __FUNCT__ 11359123193aSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqDense" 11369123193aSHong Zhang PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C) 11379123193aSHong Zhang { 11389123193aSHong Zhang PetscErrorCode ierr; 11399123193aSHong Zhang 11409123193aSHong Zhang PetscFunctionBegin; 11415a586d82SBarry Smith ierr = MatMatMultSymbolic_SeqDense_SeqDense(A,B,0.0,C);CHKERRQ(ierr); 1142d73949e8SHong Zhang (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqDense; 11439123193aSHong Zhang PetscFunctionReturn(0); 11449123193aSHong Zhang } 11459123193aSHong Zhang 11469123193aSHong Zhang #undef __FUNCT__ 11479123193aSHong Zhang #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqDense" 11489123193aSHong Zhang PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqDense(Mat A,Mat B,Mat C) 11499123193aSHong Zhang { 1150f32d5d43SBarry Smith Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data; 11519123193aSHong Zhang PetscErrorCode ierr; 1152dd6ea824SBarry Smith PetscScalar *b,*c,r1,r2,r3,r4,*b1,*b2,*b3,*b4; 1153dd6ea824SBarry Smith MatScalar *aa; 1154d0f46423SBarry Smith PetscInt cm=C->rmap->n, cn=B->cmap->n, bm=B->rmap->n, col, i,j,n,*aj, am = A->rmap->n; 1155f32d5d43SBarry Smith PetscInt am2 = 2*am, am3 = 3*am, bm4 = 4*bm,colam; 11569123193aSHong Zhang 11579123193aSHong Zhang PetscFunctionBegin; 1158f32d5d43SBarry Smith if (!cm || !cn) PetscFunctionReturn(0); 1159e32f2f54SBarry 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); 1160e32f2f54SBarry 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); 1161e32f2f54SBarry 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); 11628c778c55SBarry Smith ierr = MatDenseGetArray(B,&b);CHKERRQ(ierr); 11638c778c55SBarry Smith ierr = MatDenseGetArray(C,&c);CHKERRQ(ierr); 1164f32d5d43SBarry Smith b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm; 1165f32d5d43SBarry Smith for (col=0; col<cn-4; col += 4){ /* over columns of C */ 1166f32d5d43SBarry Smith colam = col*am; 1167f32d5d43SBarry Smith for (i=0; i<am; i++) { /* over rows of C in those columns */ 1168f32d5d43SBarry Smith r1 = r2 = r3 = r4 = 0.0; 1169f32d5d43SBarry Smith n = a->i[i+1] - a->i[i]; 1170f32d5d43SBarry Smith aj = a->j + a->i[i]; 1171f32d5d43SBarry Smith aa = a->a + a->i[i]; 1172f32d5d43SBarry Smith for (j=0; j<n; j++) { 1173f32d5d43SBarry Smith r1 += (*aa)*b1[*aj]; 1174f32d5d43SBarry Smith r2 += (*aa)*b2[*aj]; 1175f32d5d43SBarry Smith r3 += (*aa)*b3[*aj]; 1176f32d5d43SBarry Smith r4 += (*aa++)*b4[*aj++]; 11779123193aSHong Zhang } 1178f32d5d43SBarry Smith c[colam + i] = r1; 1179f32d5d43SBarry Smith c[colam + am + i] = r2; 1180f32d5d43SBarry Smith c[colam + am2 + i] = r3; 1181f32d5d43SBarry Smith c[colam + am3 + i] = r4; 1182f32d5d43SBarry Smith } 1183f32d5d43SBarry Smith b1 += bm4; 1184f32d5d43SBarry Smith b2 += bm4; 1185f32d5d43SBarry Smith b3 += bm4; 1186f32d5d43SBarry Smith b4 += bm4; 1187f32d5d43SBarry Smith } 1188f32d5d43SBarry Smith for (;col<cn; col++){ /* over extra columns of C */ 1189f32d5d43SBarry Smith for (i=0; i<am; i++) { /* over rows of C in those columns */ 1190f32d5d43SBarry Smith r1 = 0.0; 1191f32d5d43SBarry Smith n = a->i[i+1] - a->i[i]; 1192f32d5d43SBarry Smith aj = a->j + a->i[i]; 1193f32d5d43SBarry Smith aa = a->a + a->i[i]; 1194f32d5d43SBarry Smith 1195f32d5d43SBarry Smith for (j=0; j<n; j++) { 1196f32d5d43SBarry Smith r1 += (*aa++)*b1[*aj++]; 1197f32d5d43SBarry Smith } 1198f32d5d43SBarry Smith c[col*am + i] = r1; 1199f32d5d43SBarry Smith } 1200f32d5d43SBarry Smith b1 += bm; 1201f32d5d43SBarry Smith } 1202dc0b31edSSatish Balay ierr = PetscLogFlops(cn*(2.0*a->nz));CHKERRQ(ierr); 12038c778c55SBarry Smith ierr = MatDenseRestoreArray(B,&b);CHKERRQ(ierr); 12048c778c55SBarry Smith ierr = MatDenseRestoreArray(C,&c);CHKERRQ(ierr); 1205f32d5d43SBarry Smith ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1206f32d5d43SBarry Smith ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1207f32d5d43SBarry Smith PetscFunctionReturn(0); 1208f32d5d43SBarry Smith } 1209f32d5d43SBarry Smith 1210f32d5d43SBarry Smith /* 12114324174eSBarry Smith Note very similar to MatMult_SeqAIJ(), should generate both codes from same base 1212f32d5d43SBarry Smith */ 1213f32d5d43SBarry Smith #undef __FUNCT__ 1214f32d5d43SBarry Smith #define __FUNCT__ "MatMatMultNumericAdd_SeqAIJ_SeqDense" 1215f32d5d43SBarry Smith PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat A,Mat B,Mat C) 1216f32d5d43SBarry Smith { 1217f32d5d43SBarry Smith Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data; 1218f32d5d43SBarry Smith PetscErrorCode ierr; 1219dd6ea824SBarry Smith PetscScalar *b,*c,r1,r2,r3,r4,*b1,*b2,*b3,*b4; 1220dd6ea824SBarry Smith MatScalar *aa; 1221d0f46423SBarry 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; 12224324174eSBarry Smith PetscInt am2 = 2*am, am3 = 3*am, bm4 = 4*bm,colam,*ridx; 1223f32d5d43SBarry Smith 1224f32d5d43SBarry Smith PetscFunctionBegin; 1225f32d5d43SBarry Smith if (!cm || !cn) PetscFunctionReturn(0); 12268c778c55SBarry Smith ierr = MatDenseGetArray(B,&b);CHKERRQ(ierr); 12278c778c55SBarry Smith ierr = MatDenseGetArray(C,&c);CHKERRQ(ierr); 1228f32d5d43SBarry Smith b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm; 12294324174eSBarry Smith 12304324174eSBarry Smith if (a->compressedrow.use){ /* use compressed row format */ 12314324174eSBarry Smith for (col=0; col<cn-4; col += 4){ /* over columns of C */ 12324324174eSBarry Smith colam = col*am; 12334324174eSBarry Smith arm = a->compressedrow.nrows; 12344324174eSBarry Smith ii = a->compressedrow.i; 12354324174eSBarry Smith ridx = a->compressedrow.rindex; 12364324174eSBarry Smith for (i=0; i<arm; i++) { /* over rows of C in those columns */ 12374324174eSBarry Smith r1 = r2 = r3 = r4 = 0.0; 12384324174eSBarry Smith n = ii[i+1] - ii[i]; 12394324174eSBarry Smith aj = a->j + ii[i]; 12404324174eSBarry Smith aa = a->a + ii[i]; 12414324174eSBarry Smith for (j=0; j<n; j++) { 12424324174eSBarry Smith r1 += (*aa)*b1[*aj]; 12434324174eSBarry Smith r2 += (*aa)*b2[*aj]; 12444324174eSBarry Smith r3 += (*aa)*b3[*aj]; 12454324174eSBarry Smith r4 += (*aa++)*b4[*aj++]; 12464324174eSBarry Smith } 12474324174eSBarry Smith c[colam + ridx[i]] += r1; 12484324174eSBarry Smith c[colam + am + ridx[i]] += r2; 12494324174eSBarry Smith c[colam + am2 + ridx[i]] += r3; 12504324174eSBarry Smith c[colam + am3 + ridx[i]] += r4; 12514324174eSBarry Smith } 12524324174eSBarry Smith b1 += bm4; 12534324174eSBarry Smith b2 += bm4; 12544324174eSBarry Smith b3 += bm4; 12554324174eSBarry Smith b4 += bm4; 12564324174eSBarry Smith } 12574324174eSBarry Smith for (;col<cn; col++){ /* over extra columns of C */ 12584324174eSBarry Smith colam = col*am; 12594324174eSBarry Smith arm = a->compressedrow.nrows; 12604324174eSBarry Smith ii = a->compressedrow.i; 12614324174eSBarry Smith ridx = a->compressedrow.rindex; 12624324174eSBarry Smith for (i=0; i<arm; i++) { /* over rows of C in those columns */ 12634324174eSBarry Smith r1 = 0.0; 12644324174eSBarry Smith n = ii[i+1] - ii[i]; 12654324174eSBarry Smith aj = a->j + ii[i]; 12664324174eSBarry Smith aa = a->a + ii[i]; 12674324174eSBarry Smith 12684324174eSBarry Smith for (j=0; j<n; j++) { 12694324174eSBarry Smith r1 += (*aa++)*b1[*aj++]; 12704324174eSBarry Smith } 1271*a2ea699eSBarry Smith c[colam + ridx[i]] += r1; 12724324174eSBarry Smith } 12734324174eSBarry Smith b1 += bm; 12744324174eSBarry Smith } 12754324174eSBarry Smith } else { 1276f32d5d43SBarry Smith for (col=0; col<cn-4; col += 4){ /* over columns of C */ 1277f32d5d43SBarry Smith colam = col*am; 1278f32d5d43SBarry Smith for (i=0; i<am; i++) { /* over rows of C in those columns */ 1279f32d5d43SBarry Smith r1 = r2 = r3 = r4 = 0.0; 1280f32d5d43SBarry Smith n = a->i[i+1] - a->i[i]; 1281f32d5d43SBarry Smith aj = a->j + a->i[i]; 1282f32d5d43SBarry Smith aa = a->a + a->i[i]; 1283f32d5d43SBarry Smith for (j=0; j<n; j++) { 1284f32d5d43SBarry Smith r1 += (*aa)*b1[*aj]; 1285f32d5d43SBarry Smith r2 += (*aa)*b2[*aj]; 1286f32d5d43SBarry Smith r3 += (*aa)*b3[*aj]; 1287f32d5d43SBarry Smith r4 += (*aa++)*b4[*aj++]; 1288f32d5d43SBarry Smith } 1289f32d5d43SBarry Smith c[colam + i] += r1; 1290f32d5d43SBarry Smith c[colam + am + i] += r2; 1291f32d5d43SBarry Smith c[colam + am2 + i] += r3; 1292f32d5d43SBarry Smith c[colam + am3 + i] += r4; 1293f32d5d43SBarry Smith } 1294f32d5d43SBarry Smith b1 += bm4; 1295f32d5d43SBarry Smith b2 += bm4; 1296f32d5d43SBarry Smith b3 += bm4; 1297f32d5d43SBarry Smith b4 += bm4; 1298f32d5d43SBarry Smith } 1299f32d5d43SBarry Smith for (;col<cn; col++){ /* over extra columns of C */ 1300*a2ea699eSBarry Smith colam = col*am; 1301f32d5d43SBarry Smith for (i=0; i<am; i++) { /* over rows of C in those columns */ 1302f32d5d43SBarry Smith r1 = 0.0; 1303f32d5d43SBarry Smith n = a->i[i+1] - a->i[i]; 1304f32d5d43SBarry Smith aj = a->j + a->i[i]; 1305f32d5d43SBarry Smith aa = a->a + a->i[i]; 1306f32d5d43SBarry Smith 1307f32d5d43SBarry Smith for (j=0; j<n; j++) { 1308f32d5d43SBarry Smith r1 += (*aa++)*b1[*aj++]; 1309f32d5d43SBarry Smith } 1310*a2ea699eSBarry Smith c[colam + i] += r1; 1311f32d5d43SBarry Smith } 1312f32d5d43SBarry Smith b1 += bm; 1313f32d5d43SBarry Smith } 13144324174eSBarry Smith } 1315dc0b31edSSatish Balay ierr = PetscLogFlops(cn*2.0*a->nz);CHKERRQ(ierr); 13168c778c55SBarry Smith ierr = MatDenseRestoreArray(B,&b);CHKERRQ(ierr); 13178c778c55SBarry Smith ierr = MatDenseRestoreArray(C,&c);CHKERRQ(ierr); 13189123193aSHong Zhang PetscFunctionReturn(0); 13199123193aSHong Zhang } 1320b1683b59SHong Zhang 1321b1683b59SHong Zhang #undef __FUNCT__ 1322b9af6bddSHong Zhang #define __FUNCT__ "MatTransColoringApplySpToDen_SeqAIJ" 1323b9af6bddSHong Zhang PetscErrorCode MatTransColoringApplySpToDen_SeqAIJ(MatTransposeColoring coloring,Mat B,Mat Btdense) 1324c8db22f6SHong Zhang { 1325c8db22f6SHong Zhang PetscErrorCode ierr; 13262f5f1f90SHong Zhang Mat_SeqAIJ *b = (Mat_SeqAIJ*)B->data; 13272f5f1f90SHong Zhang Mat_SeqDense *btdense = (Mat_SeqDense*)Btdense->data; 13282f5f1f90SHong Zhang PetscInt *bi=b->i,*bj=b->j; 13292f5f1f90SHong Zhang PetscInt m=Btdense->rmap->n,n=Btdense->cmap->n,j,k,l,col,anz,*btcol,brow,ncolumns; 13302f5f1f90SHong Zhang MatScalar *btval,*btval_den,*ba=b->a; 13312f5f1f90SHong Zhang PetscInt *columns=coloring->columns,*colorforcol=coloring->colorforcol,ncolors=coloring->ncolors; 1332c8db22f6SHong Zhang 1333c8db22f6SHong Zhang PetscFunctionBegin; 13342f5f1f90SHong Zhang btval_den=btdense->v; 13352f5f1f90SHong Zhang ierr = PetscMemzero(btval_den,(m*n)*sizeof(MatScalar));CHKERRQ(ierr); 13362f5f1f90SHong Zhang for (k=0; k<ncolors; k++) { 13372f5f1f90SHong Zhang ncolumns = coloring->ncolumns[k]; 13382f5f1f90SHong Zhang for (l=0; l<ncolumns; l++) { /* insert a row of B to a column of Btdense */ 1339d2853540SHong Zhang col = *(columns + colorforcol[k] + l); 13402f5f1f90SHong Zhang btcol = bj + bi[col]; 13412f5f1f90SHong Zhang btval = ba + bi[col]; 13422f5f1f90SHong Zhang anz = bi[col+1] - bi[col]; 1343d2853540SHong Zhang for (j=0; j<anz; j++){ 13442f5f1f90SHong Zhang brow = btcol[j]; 13452f5f1f90SHong Zhang btval_den[brow] = btval[j]; 1346c8db22f6SHong Zhang } 1347c8db22f6SHong Zhang } 13482f5f1f90SHong Zhang btval_den += m; 1349c8db22f6SHong Zhang } 1350c8db22f6SHong Zhang PetscFunctionReturn(0); 1351c8db22f6SHong Zhang } 1352c8db22f6SHong Zhang 1353c8db22f6SHong Zhang #undef __FUNCT__ 1354b9af6bddSHong Zhang #define __FUNCT__ "MatTransColoringApplyDenToSp_SeqAIJ" 1355b9af6bddSHong Zhang PetscErrorCode MatTransColoringApplyDenToSp_SeqAIJ(MatTransposeColoring matcoloring,Mat Cden,Mat Csp) 13568972f759SHong Zhang { 13578972f759SHong Zhang PetscErrorCode ierr; 1358b2d2b04fSHong Zhang Mat_SeqAIJ *csp = (Mat_SeqAIJ*)Csp->data; 13592f5f1f90SHong Zhang PetscInt k,l,*row,*idx,m,ncolors=matcoloring->ncolors,nrows; 1360b2d2b04fSHong Zhang PetscScalar *ca_den,*cp_den,*ca=csp->a; 13612f5f1f90SHong Zhang PetscInt *rows=matcoloring->rows,*spidx=matcoloring->columnsforspidx,*colorforrow=matcoloring->colorforrow; 13628972f759SHong Zhang 13638972f759SHong Zhang PetscFunctionBegin; 13648972f759SHong Zhang ierr = MatGetLocalSize(Csp,&m,PETSC_NULL);CHKERRQ(ierr); 1365a3fe58edSHong Zhang ierr = MatDenseGetArray(Cden,&ca_den);CHKERRQ(ierr); 1366b2d2b04fSHong Zhang cp_den = ca_den; 13672f5f1f90SHong Zhang for (k=0; k<ncolors; k++) { 13682f5f1f90SHong Zhang nrows = matcoloring->nrows[k]; 13692f5f1f90SHong Zhang row = rows + colorforrow[k]; 13702f5f1f90SHong Zhang idx = spidx + colorforrow[k]; 13712f5f1f90SHong Zhang for (l=0; l<nrows; l++){ 13722f5f1f90SHong Zhang ca[idx[l]] = cp_den[row[l]]; 1373b2d2b04fSHong Zhang } 1374b2d2b04fSHong Zhang cp_den += m; 1375b2d2b04fSHong Zhang } 1376a3fe58edSHong Zhang ierr = MatDenseRestoreArray(Cden,&ca_den);CHKERRQ(ierr); 13778972f759SHong Zhang PetscFunctionReturn(0); 13788972f759SHong Zhang } 13798972f759SHong Zhang 1380e99be685SHong Zhang /* 1381e99be685SHong Zhang MatGetColumnIJ_SeqAIJ_Color() and MatRestoreColumnIJ_SeqAIJ_Color() are customized from 1382e99be685SHong Zhang MatGetColumnIJ_SeqAIJ() and MatRestoreColumnIJ_SeqAIJ() by adding an output 1383e99be685SHong Zhang spidx[], index of a->j, to be used for setting 'columnsforspidx' in MatTransposeColoringCreate_SeqAIJ(). 1384e99be685SHong Zhang */ 1385e99be685SHong Zhang #undef __FUNCT__ 1386e99be685SHong Zhang #define __FUNCT__ "MatGetColumnIJ_SeqAIJ_Color" 13871a83f524SJed Brown PetscErrorCode MatGetColumnIJ_SeqAIJ_Color(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool inodecompressed,PetscInt *nn,const PetscInt *ia[],const PetscInt *ja[],PetscInt *spidx[],PetscBool *done) 1388e99be685SHong Zhang { 1389e99be685SHong Zhang Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 1390e99be685SHong Zhang PetscErrorCode ierr; 1391e99be685SHong Zhang PetscInt i,*collengths,*cia,*cja,n = A->cmap->n,m = A->rmap->n; 1392e99be685SHong Zhang PetscInt nz = a->i[m],row,*jj,mr,col; 1393e99be685SHong Zhang PetscInt *cspidx; 1394e99be685SHong Zhang 1395e99be685SHong Zhang PetscFunctionBegin; 1396e99be685SHong Zhang *nn = n; 1397e99be685SHong Zhang if (!ia) PetscFunctionReturn(0); 1398e99be685SHong Zhang if (symmetric) { 1399e99be685SHong Zhang SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"MatGetColumnIJ_SeqAIJ_Color() not supported for the case symmetric"); 14001a83f524SJed Brown ierr = MatToSymmetricIJ_SeqAIJ(A->rmap->n,a->i,a->j,0,oshift,(PetscInt**)ia,(PetscInt**)ja);CHKERRQ(ierr); 1401e99be685SHong Zhang } else { 1402e99be685SHong Zhang ierr = PetscMalloc((n+1)*sizeof(PetscInt),&collengths);CHKERRQ(ierr); 1403e99be685SHong Zhang ierr = PetscMemzero(collengths,n*sizeof(PetscInt));CHKERRQ(ierr); 1404e99be685SHong Zhang ierr = PetscMalloc((n+1)*sizeof(PetscInt),&cia);CHKERRQ(ierr); 1405e99be685SHong Zhang ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&cja);CHKERRQ(ierr); 1406e99be685SHong Zhang ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&cspidx);CHKERRQ(ierr); 1407e99be685SHong Zhang jj = a->j; 1408e99be685SHong Zhang for (i=0; i<nz; i++) { 1409e99be685SHong Zhang collengths[jj[i]]++; 1410e99be685SHong Zhang } 1411e99be685SHong Zhang cia[0] = oshift; 1412e99be685SHong Zhang for (i=0; i<n; i++) { 1413e99be685SHong Zhang cia[i+1] = cia[i] + collengths[i]; 1414e99be685SHong Zhang } 1415e99be685SHong Zhang ierr = PetscMemzero(collengths,n*sizeof(PetscInt));CHKERRQ(ierr); 1416e99be685SHong Zhang jj = a->j; 1417e99be685SHong Zhang for (row=0; row<m; row++) { 1418e99be685SHong Zhang mr = a->i[row+1] - a->i[row]; 1419e99be685SHong Zhang for (i=0; i<mr; i++) { 1420e99be685SHong Zhang col = *jj++; 1421e99be685SHong Zhang cspidx[cia[col] + collengths[col] - oshift] = a->i[row] + i; /* index of a->j */ 1422e99be685SHong Zhang cja[cia[col] + collengths[col]++ - oshift] = row + oshift; 1423e99be685SHong Zhang } 1424e99be685SHong Zhang } 1425e99be685SHong Zhang ierr = PetscFree(collengths);CHKERRQ(ierr); 1426e99be685SHong Zhang *ia = cia; *ja = cja; 1427e99be685SHong Zhang *spidx = cspidx; 1428e99be685SHong Zhang } 1429e99be685SHong Zhang PetscFunctionReturn(0); 1430e99be685SHong Zhang } 1431e99be685SHong Zhang 1432e99be685SHong Zhang #undef __FUNCT__ 1433e99be685SHong Zhang #define __FUNCT__ "MatRestoreColumnIJ_SeqAIJ_Color" 14341a83f524SJed Brown PetscErrorCode MatRestoreColumnIJ_SeqAIJ_Color(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool inodecompressed,PetscInt *n,const PetscInt *ia[],const PetscInt *ja[],PetscInt *spidx[],PetscBool *done) 1435e99be685SHong Zhang { 1436e99be685SHong Zhang PetscErrorCode ierr; 1437e99be685SHong Zhang 1438e99be685SHong Zhang PetscFunctionBegin; 1439e99be685SHong Zhang if (!ia) PetscFunctionReturn(0); 1440e99be685SHong Zhang 1441e99be685SHong Zhang ierr = PetscFree(*ia);CHKERRQ(ierr); 1442e99be685SHong Zhang ierr = PetscFree(*ja);CHKERRQ(ierr); 1443e99be685SHong Zhang ierr = PetscFree(*spidx);CHKERRQ(ierr); 1444e99be685SHong Zhang PetscFunctionReturn(0); 1445e99be685SHong Zhang } 1446e99be685SHong Zhang 14478972f759SHong Zhang #undef __FUNCT__ 1448b9af6bddSHong Zhang #define __FUNCT__ "MatTransposeColoringCreate_SeqAIJ" 1449b9af6bddSHong Zhang PetscErrorCode MatTransposeColoringCreate_SeqAIJ(Mat mat,ISColoring iscoloring,MatTransposeColoring c) 1450b1683b59SHong Zhang { 1451b1683b59SHong Zhang PetscErrorCode ierr; 14521a83f524SJed Brown PetscInt i,n,nrows,N,j,k,m,ncols,col,cm; 14531a83f524SJed Brown const PetscInt *is,*ci,*cj,*row_idx; 1454b28d80bdSHong Zhang PetscInt nis = iscoloring->n,*rowhit,bs = 1; 1455b1683b59SHong Zhang IS *isa; 1456b1683b59SHong Zhang PetscBool flg1,flg2; 1457b28d80bdSHong Zhang Mat_SeqAIJ *csp = (Mat_SeqAIJ*)mat->data; 1458e99be685SHong Zhang PetscInt *colorforrow,*rows,*rows_i,*columnsforspidx,*columnsforspidx_i,*idxhit,*spidx; 1459d2853540SHong Zhang PetscInt *colorforcol,*columns,*columns_i; 1460b1683b59SHong Zhang 1461b1683b59SHong Zhang PetscFunctionBegin; 1462b1683b59SHong Zhang ierr = ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);CHKERRQ(ierr); 1463e99be685SHong Zhang 1464b1683b59SHong Zhang /* this is ugly way to get blocksize but cannot call MatGetBlockSize() because AIJ can have bs > 1 */ 1465251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)mat,MATSEQBAIJ,&flg1);CHKERRQ(ierr); 1466251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)mat,MATMPIBAIJ,&flg2);CHKERRQ(ierr); 1467b1683b59SHong Zhang if (flg1 || flg2) { 1468b1683b59SHong Zhang ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr); 1469b1683b59SHong Zhang } 1470b1683b59SHong Zhang 1471b1683b59SHong Zhang N = mat->cmap->N/bs; 1472b1683b59SHong Zhang c->M = mat->rmap->N/bs; /* set total rows, columns and local rows */ 1473b1683b59SHong Zhang c->N = mat->cmap->N/bs; 1474b1683b59SHong Zhang c->m = mat->rmap->N/bs; 1475b1683b59SHong Zhang c->rstart = 0; 1476b1683b59SHong Zhang 1477b1683b59SHong Zhang c->ncolors = nis; 1478b1683b59SHong Zhang ierr = PetscMalloc(nis*sizeof(PetscInt),&c->ncolumns);CHKERRQ(ierr); 1479b28d80bdSHong Zhang ierr = PetscMalloc(nis*sizeof(PetscInt),&c->nrows);CHKERRQ(ierr); 1480d2853540SHong Zhang ierr = PetscMalloc2(csp->nz+1,PetscInt,&rows,csp->nz+1,PetscInt,&columnsforspidx);CHKERRQ(ierr); 1481d2853540SHong Zhang ierr = PetscMalloc((nis+1)*sizeof(PetscInt),&colorforrow);CHKERRQ(ierr); 1482d2853540SHong Zhang colorforrow[0] = 0; 1483d2853540SHong Zhang rows_i = rows; 1484d2853540SHong Zhang columnsforspidx_i = columnsforspidx; 1485b1683b59SHong Zhang 1486d2853540SHong Zhang ierr = PetscMalloc((nis+1)*sizeof(PetscInt),&colorforcol);CHKERRQ(ierr); 1487d2853540SHong Zhang ierr = PetscMalloc((N+1)*sizeof(PetscInt),&columns);CHKERRQ(ierr); 1488d2853540SHong Zhang colorforcol[0] = 0; 1489d2853540SHong Zhang columns_i = columns; 1490d2853540SHong Zhang 1491f2691f8dSHong Zhang ierr = MatGetColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,PETSC_NULL);CHKERRQ(ierr); /* column-wise storage of mat */ 1492b1683b59SHong Zhang 1493b28d80bdSHong Zhang cm = c->m; 1494b28d80bdSHong Zhang ierr = PetscMalloc((cm+1)*sizeof(PetscInt),&rowhit);CHKERRQ(ierr); 1495e99be685SHong Zhang ierr = PetscMalloc((cm+1)*sizeof(PetscInt),&idxhit);CHKERRQ(ierr); 1496b1683b59SHong Zhang for (i=0; i<nis; i++) { 1497b1683b59SHong Zhang ierr = ISGetLocalSize(isa[i],&n);CHKERRQ(ierr); 1498b1683b59SHong Zhang ierr = ISGetIndices(isa[i],&is);CHKERRQ(ierr); 1499b1683b59SHong Zhang c->ncolumns[i] = n; 1500b1683b59SHong Zhang if (n) { 1501d2853540SHong Zhang ierr = PetscMemcpy(columns_i,is,n*sizeof(PetscInt));CHKERRQ(ierr); 1502b1683b59SHong Zhang } 1503d2853540SHong Zhang colorforcol[i+1] = colorforcol[i] + n; 1504d2853540SHong Zhang columns_i += n; 1505b1683b59SHong Zhang 1506b1683b59SHong Zhang /* fast, crude version requires O(N*N) work */ 1507e8354b3eSHong Zhang ierr = PetscMemzero(rowhit,cm*sizeof(PetscInt));CHKERRQ(ierr); 1508e99be685SHong Zhang 1509b1683b59SHong Zhang /* loop over columns*/ 1510b1683b59SHong Zhang for (j=0; j<n; j++) { 1511b1683b59SHong Zhang col = is[j]; 1512d2853540SHong Zhang row_idx = cj + ci[col]; 1513b1683b59SHong Zhang m = ci[col+1] - ci[col]; 1514b1683b59SHong Zhang /* loop over columns marking them in rowhit */ 1515b1683b59SHong Zhang for (k=0; k<m; k++) { 1516e99be685SHong Zhang idxhit[*row_idx] = spidx[ci[col] + k]; 1517d2853540SHong Zhang rowhit[*row_idx++] = col + 1; 1518b1683b59SHong Zhang } 1519b1683b59SHong Zhang } 1520b1683b59SHong Zhang /* count the number of hits */ 1521b1683b59SHong Zhang nrows = 0; 1522e8354b3eSHong Zhang for (j=0; j<cm; j++) { 1523b1683b59SHong Zhang if (rowhit[j]) nrows++; 1524b1683b59SHong Zhang } 1525b1683b59SHong Zhang c->nrows[i] = nrows; 1526d2853540SHong Zhang colorforrow[i+1] = colorforrow[i] + nrows; 1527d2853540SHong Zhang 1528b1683b59SHong Zhang nrows = 0; 1529e8354b3eSHong Zhang for (j=0; j<cm; j++) { 1530b1683b59SHong Zhang if (rowhit[j]) { 1531d2853540SHong Zhang rows_i[nrows] = j; 1532e99be685SHong Zhang columnsforspidx_i[nrows] = idxhit[j]; 1533b1683b59SHong Zhang nrows++; 1534b1683b59SHong Zhang } 1535b1683b59SHong Zhang } 1536b1683b59SHong Zhang ierr = ISRestoreIndices(isa[i],&is);CHKERRQ(ierr); 1537d2853540SHong Zhang rows_i += nrows; columnsforspidx_i += nrows; 1538b1683b59SHong Zhang } 1539f2691f8dSHong Zhang ierr = MatRestoreColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,PETSC_NULL);CHKERRQ(ierr); 1540b28d80bdSHong Zhang ierr = PetscFree(rowhit);CHKERRQ(ierr); 1541b1683b59SHong Zhang ierr = ISColoringRestoreIS(iscoloring,&isa);CHKERRQ(ierr); 1542d2853540SHong Zhang #if defined(PETSC_USE_DEBUG) 1543d2853540SHong Zhang if (csp->nz != colorforrow[nis]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"csp->nz %d != colorforrow[nis] %d",csp->nz,colorforrow[nis]); 1544d2853540SHong Zhang #endif 1545b28d80bdSHong Zhang 1546d2853540SHong Zhang c->colorforrow = colorforrow; 1547d2853540SHong Zhang c->rows = rows; 1548d2853540SHong Zhang c->columnsforspidx = columnsforspidx; 1549d2853540SHong Zhang c->colorforcol = colorforcol; 1550d2853540SHong Zhang c->columns = columns; 1551f94ccd6cSHong Zhang ierr = PetscFree(idxhit);CHKERRQ(ierr); 1552b1683b59SHong Zhang PetscFunctionReturn(0); 1553b1683b59SHong Zhang } 1554