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*/ 12e9e4536cSHong Zhang /* 13e9e4536cSHong Zhang #define DEBUG_MATMATMULT 14e9e4536cSHong Zhang */ 15ec01deb9SMatthew Knepley EXTERN_C_BEGIN 166284ec50SHong Zhang #undef __FUNCT__ 175c66b693SKris Buschelman #define __FUNCT__ "MatMatMult_SeqAIJ_SeqAIJ" 1838baddfdSBarry Smith PetscErrorCode MatMatMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 1938baddfdSBarry Smith { 20dfbe8321SBarry Smith PetscErrorCode ierr; 218a07c6f1SJed Brown PetscBool scalable=PETSC_FALSE,scalable_fast=PETSC_FALSE,heap = PETSC_FALSE,btheap = PETSC_FALSE; 225c66b693SKris Buschelman 235c66b693SKris Buschelman PetscFunctionBegin; 2426be0446SHong Zhang if (scall == MAT_INITIAL_MATRIX){ 25152983bfSHong Zhang ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr); 26152983bfSHong Zhang ierr = PetscOptionsBool("-matmatmult_scalable","Use a scalable but slower C=A*B","",scalable,&scalable,PETSC_NULL);CHKERRQ(ierr); 273c50cad2SHong Zhang ierr = PetscOptionsBool("-matmatmult_scalable_fast","Use a scalable but slower C=A*B","",scalable_fast,&scalable_fast,PETSC_NULL);CHKERRQ(ierr); 280ced3a2bSJed Brown ierr = PetscOptionsBool("-matmatmult_heap","Use heap implementation of symbolic factorization C=A*B","",heap,&heap,PETSC_NULL);CHKERRQ(ierr); 298a07c6f1SJed Brown ierr = PetscOptionsBool("-matmatmult_btheap","Use btheap implementation of symbolic factorization C=A*B","",btheap,&btheap,PETSC_NULL);CHKERRQ(ierr); 30d8bbc50fSBarry Smith ierr = PetscOptionsEnd();CHKERRQ(ierr); 31d8bbc50fSBarry Smith ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr); 323c50cad2SHong Zhang if (scalable_fast){ 333c50cad2SHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast(A,B,fill,C);CHKERRQ(ierr); 34aacf854cSBarry Smith } else if (scalable){ 35aacf854cSBarry Smith ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(A,B,fill,C);CHKERRQ(ierr); 360ced3a2bSJed Brown } else if (heap) { 370ced3a2bSJed Brown ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_Heap(A,B,fill,C);CHKERRQ(ierr); 388a07c6f1SJed Brown } else if (btheap) { 398a07c6f1SJed Brown ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_BTHeap(A,B,fill,C);CHKERRQ(ierr); 4025023028SHong Zhang } else { 4126be0446SHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr); 4225023028SHong Zhang } 43598bc09dSHong Zhang ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr); 4426be0446SHong Zhang } 45598bc09dSHong Zhang ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr); 4601e47db4SHong Zhang ierr = (*(*C)->ops->matmultnumeric)(A,B,*C);CHKERRQ(ierr); 47598bc09dSHong Zhang ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr); 485c66b693SKris Buschelman PetscFunctionReturn(0); 495c66b693SKris Buschelman } 50ec01deb9SMatthew Knepley EXTERN_C_END 511c24bd37SHong Zhang 52e9e4536cSHong Zhang #undef __FUNCT__ 53b561aa9dSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ" 54b561aa9dSHong Zhang PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C) 55b561aa9dSHong Zhang { 56b561aa9dSHong Zhang PetscErrorCode ierr; 57b561aa9dSHong Zhang Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c; 58971236abSHong Zhang PetscInt *ai=a->i,*bi=b->i,*ci,*cj; 59c1ba5575SJed Brown PetscInt am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N; 60b561aa9dSHong Zhang PetscReal afill; 61fb908031SHong Zhang PetscInt i,j,anzi,brow,bnzj,cnzi,*bj,*aj,nlnk_max,*lnk,ndouble=0; 62fb908031SHong Zhang PetscBT lnkbt; 63fb908031SHong Zhang PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 64b561aa9dSHong Zhang 65b561aa9dSHong Zhang PetscFunctionBegin; 66bd958071SHong Zhang /* Get ci and cj */ 67bd958071SHong Zhang /*---------------*/ 68fb908031SHong Zhang /* Allocate ci array, arrays for fill computation and */ 69fb908031SHong Zhang /* free space for accumulating nonzero column info */ 70fb908031SHong Zhang ierr = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr); 71fb908031SHong Zhang ci[0] = 0; 72fb908031SHong Zhang 73fb908031SHong Zhang /* create and initialize a linked list */ 74fb908031SHong Zhang nlnk_max = a->rmax*b->rmax; 75fb908031SHong Zhang if (!nlnk_max || nlnk_max > bn) nlnk_max = bn; 76fb908031SHong Zhang ierr = PetscLLCondensedCreate(nlnk_max,bn,&lnk,&lnkbt);CHKERRQ(ierr); 77fb908031SHong Zhang 78fb908031SHong Zhang /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */ 79fb908031SHong Zhang ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr); 80fb908031SHong Zhang current_space = free_space; 81fb908031SHong Zhang 82fb908031SHong Zhang /* Determine ci and cj */ 83fb908031SHong Zhang for (i=0; i<am; i++) { 84fb908031SHong Zhang anzi = ai[i+1] - ai[i]; 85fb908031SHong Zhang aj = a->j + ai[i]; 86fb908031SHong Zhang for (j=0; j<anzi; j++){ 87fb908031SHong Zhang brow = aj[j]; 88fb908031SHong Zhang bnzj = bi[brow+1] - bi[brow]; 89fb908031SHong Zhang bj = b->j + bi[brow]; 90fb908031SHong Zhang /* add non-zero cols of B into the sorted linked list lnk */ 91fb908031SHong Zhang ierr = PetscLLCondensedAddSorted(bnzj,bj,lnk,lnkbt);CHKERRQ(ierr); 92fb908031SHong Zhang } 93fb908031SHong Zhang cnzi = lnk[0]; 94fb908031SHong Zhang 95fb908031SHong Zhang /* If free space is not available, make more free space */ 96fb908031SHong Zhang /* Double the amount of total space in the list */ 97fb908031SHong Zhang if (current_space->local_remaining<cnzi) { 98fb908031SHong Zhang ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 99fb908031SHong Zhang ndouble++; 100fb908031SHong Zhang } 101fb908031SHong Zhang 102fb908031SHong Zhang /* Copy data into free space, then initialize lnk */ 103fb908031SHong Zhang ierr = PetscLLCondensedClean(bn,cnzi,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 104fb908031SHong Zhang current_space->array += cnzi; 105fb908031SHong Zhang current_space->local_used += cnzi; 106fb908031SHong Zhang current_space->local_remaining -= cnzi; 107fb908031SHong Zhang ci[i+1] = ci[i] + cnzi; 108fb908031SHong Zhang } 109fb908031SHong Zhang 110fb908031SHong Zhang /* Column indices are in the list of free space */ 111fb908031SHong Zhang /* Allocate space for cj, initialize cj, and */ 112fb908031SHong Zhang /* destroy list of free space and other temporary array(s) */ 113fb908031SHong Zhang ierr = PetscMalloc((ci[am]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr); 114fb908031SHong Zhang ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr); 115fb908031SHong Zhang ierr = PetscLLCondensedDestroy(lnk,lnkbt);CHKERRQ(ierr); 116b561aa9dSHong Zhang 11726be0446SHong Zhang /* put together the new symbolic matrix */ 118aa1aec99SHong Zhang ierr = MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,am,bn,ci,cj,PETSC_NULL,C);CHKERRQ(ierr); 119a2f3521dSMark F. Adams (*C)->rmap->bs = A->rmap->bs; 120a2f3521dSMark F. Adams (*C)->cmap->bs = B->cmap->bs; 12158c24d83SHong Zhang 12258c24d83SHong Zhang /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 12358c24d83SHong Zhang /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */ 12458c24d83SHong Zhang c = (Mat_SeqAIJ *)((*C)->data); 125aa1aec99SHong Zhang c->free_a = PETSC_FALSE; 126e6b907acSBarry Smith c->free_ij = PETSC_TRUE; 12758c24d83SHong Zhang c->nonew = 0; 128002e8affSHong Zhang (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ; /* fast, needs non-scalable O(bn) array 'abdense' */ 1290b7e3e3dSHong Zhang 1307212da7cSHong Zhang /* set MatInfo */ 1317212da7cSHong Zhang afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5; 132f2b054eeSHong Zhang if (afill < 1.0) afill = 1.0; 1337212da7cSHong Zhang c->maxnz = ci[am]; 1347212da7cSHong Zhang c->nz = ci[am]; 135fb908031SHong Zhang (*C)->info.mallocs = ndouble; 1367212da7cSHong Zhang (*C)->info.fill_ratio_given = fill; 1377212da7cSHong Zhang (*C)->info.fill_ratio_needed = afill; 1387212da7cSHong Zhang 1397212da7cSHong Zhang #if defined(PETSC_USE_INFO) 1407212da7cSHong Zhang if (ci[am]) { 141fb908031SHong Zhang ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",ndouble,fill,afill);CHKERRQ(ierr); 142f2b054eeSHong Zhang ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr); 143f2b054eeSHong Zhang } else { 144f2b054eeSHong Zhang ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr); 145be0fcf8dSHong Zhang } 146f2b054eeSHong Zhang #endif 14758c24d83SHong Zhang PetscFunctionReturn(0); 14858c24d83SHong Zhang } 149d50806bdSBarry Smith 15026be0446SHong Zhang #undef __FUNCT__ 15126be0446SHong Zhang #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqAIJ" 152dfbe8321SBarry Smith PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C) 153d50806bdSBarry Smith { 154dfbe8321SBarry Smith PetscErrorCode ierr; 155d13dce4bSSatish Balay PetscLogDouble flops=0.0; 156d50806bdSBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 157d50806bdSBarry Smith Mat_SeqAIJ *b = (Mat_SeqAIJ *)B->data; 158d50806bdSBarry Smith Mat_SeqAIJ *c = (Mat_SeqAIJ *)C->data; 15938baddfdSBarry Smith PetscInt *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j; 160c58ca2e3SHong Zhang PetscInt am=A->rmap->n,cm=C->rmap->n; 161519eb980SHong Zhang PetscInt i,j,k,anzi,bnzi,cnzi,brow; 162aa1aec99SHong Zhang PetscScalar *aa=a->a,*ba=b->a,*baj,*ca,valtmp; 163aa1aec99SHong Zhang PetscScalar *ab_dense; 164d50806bdSBarry Smith 165d50806bdSBarry Smith PetscFunctionBegin; 166aa1aec99SHong Zhang /* printf("MatMatMultNumeric_SeqAIJ_SeqAIJ...ca %p\n",c->a); */ 167aa1aec99SHong Zhang if (!c->a){ /* first call of MatMatMultNumeric_SeqAIJ_SeqAIJ, allocate ca and matmult_abdense */ 168aa1aec99SHong Zhang ierr = PetscMalloc((ci[cm]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr); 169aa1aec99SHong Zhang c->a = ca; 170aa1aec99SHong Zhang c->free_a = PETSC_TRUE; 171aa1aec99SHong Zhang 172aa1aec99SHong Zhang ierr = PetscMalloc(B->cmap->N*sizeof(PetscScalar),&ab_dense);CHKERRQ(ierr); 173aa1aec99SHong Zhang ierr = PetscMemzero(ab_dense,B->cmap->N*sizeof(PetscScalar));CHKERRQ(ierr); 174aa1aec99SHong Zhang c->matmult_abdense = ab_dense; 175aa1aec99SHong Zhang } else { 176aa1aec99SHong Zhang ca = c->a; 177aa1aec99SHong Zhang ab_dense = c->matmult_abdense; 178aa1aec99SHong Zhang } 179aa1aec99SHong Zhang 180c124e916SHong Zhang /* clean old values in C */ 181c124e916SHong Zhang ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr); 182d50806bdSBarry Smith /* Traverse A row-wise. */ 183d50806bdSBarry Smith /* Build the ith row in C by summing over nonzero columns in A, */ 184d50806bdSBarry Smith /* the rows of B corresponding to nonzeros of A. */ 185d50806bdSBarry Smith for (i=0; i<am; i++) { 186d50806bdSBarry Smith anzi = ai[i+1] - ai[i]; 187d50806bdSBarry Smith for (j=0; j<anzi; j++) { 188519eb980SHong Zhang brow = aj[j]; 189d50806bdSBarry Smith bnzi = bi[brow+1] - bi[brow]; 190d50806bdSBarry Smith bjj = bj + bi[brow]; 191d50806bdSBarry Smith baj = ba + bi[brow]; 192519eb980SHong Zhang /* perform dense axpy */ 19336ec6d2dSHong Zhang valtmp = aa[j]; 194519eb980SHong Zhang for (k=0; k<bnzi; k++) { 19536ec6d2dSHong Zhang ab_dense[bjj[k]] += valtmp*baj[k]; 196519eb980SHong Zhang } 197519eb980SHong Zhang flops += 2*bnzi; 198519eb980SHong Zhang } 199c58ca2e3SHong Zhang aj += anzi; aa += anzi; 200c58ca2e3SHong Zhang 201c58ca2e3SHong Zhang cnzi = ci[i+1] - ci[i]; 202519eb980SHong Zhang for (k=0; k<cnzi; k++) { 203519eb980SHong Zhang ca[k] += ab_dense[cj[k]]; 204519eb980SHong Zhang ab_dense[cj[k]] = 0.0; /* zero ab_dense */ 205519eb980SHong Zhang } 206519eb980SHong Zhang flops += cnzi; 207519eb980SHong Zhang cj += cnzi; ca += cnzi; 208519eb980SHong Zhang } 209c58ca2e3SHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 210c58ca2e3SHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 211c58ca2e3SHong Zhang ierr = PetscLogFlops(flops);CHKERRQ(ierr); 212c58ca2e3SHong Zhang PetscFunctionReturn(0); 213c58ca2e3SHong Zhang } 214c58ca2e3SHong Zhang 215c58ca2e3SHong Zhang #undef __FUNCT__ 21625023028SHong Zhang #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable" 21725023028SHong Zhang PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(Mat A,Mat B,Mat C) 218c58ca2e3SHong Zhang { 219c58ca2e3SHong Zhang PetscErrorCode ierr; 220c58ca2e3SHong Zhang PetscLogDouble flops=0.0; 221c58ca2e3SHong Zhang Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 222c58ca2e3SHong Zhang Mat_SeqAIJ *b = (Mat_SeqAIJ *)B->data; 223c58ca2e3SHong Zhang Mat_SeqAIJ *c = (Mat_SeqAIJ *)C->data; 224c58ca2e3SHong Zhang PetscInt *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j; 225c58ca2e3SHong Zhang PetscInt am=A->rmap->N,cm=C->rmap->N; 226c58ca2e3SHong Zhang PetscInt i,j,k,anzi,bnzi,cnzi,brow; 22736ec6d2dSHong Zhang PetscScalar *aa=a->a,*ba=b->a,*baj,*ca=c->a,valtmp; 228c58ca2e3SHong Zhang PetscInt nextb; 229c58ca2e3SHong Zhang 230c58ca2e3SHong Zhang PetscFunctionBegin; 231e9e4536cSHong Zhang #if defined(DEBUG_MATMATMULT) 232971236abSHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF,"MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable...\n");CHKERRQ(ierr); 233e9e4536cSHong Zhang #endif 234c58ca2e3SHong Zhang /* clean old values in C */ 235c58ca2e3SHong Zhang ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr); 236c58ca2e3SHong Zhang /* Traverse A row-wise. */ 237c58ca2e3SHong Zhang /* Build the ith row in C by summing over nonzero columns in A, */ 238c58ca2e3SHong Zhang /* the rows of B corresponding to nonzeros of A. */ 239519eb980SHong Zhang for (i=0;i<am;i++) { 240519eb980SHong Zhang anzi = ai[i+1] - ai[i]; 241519eb980SHong Zhang cnzi = ci[i+1] - ci[i]; 242519eb980SHong Zhang for (j=0;j<anzi;j++) { 243519eb980SHong Zhang brow = aj[j]; 244519eb980SHong Zhang bnzi = bi[brow+1] - bi[brow]; 245519eb980SHong Zhang bjj = bj + bi[brow]; 246519eb980SHong Zhang baj = ba + bi[brow]; 247519eb980SHong Zhang /* perform sparse axpy */ 24836ec6d2dSHong Zhang valtmp = aa[j]; 249c124e916SHong Zhang nextb = 0; 250c124e916SHong Zhang for (k=0; nextb<bnzi; k++) { 251c124e916SHong Zhang if (cj[k] == bjj[nextb]){ /* ccol == bcol */ 25236ec6d2dSHong Zhang ca[k] += valtmp*baj[nextb++]; 253c124e916SHong Zhang } 254d50806bdSBarry Smith } 255d50806bdSBarry Smith flops += 2*bnzi; 256d50806bdSBarry Smith } 257519eb980SHong Zhang aj += anzi; aa += anzi; 258519eb980SHong Zhang cj += cnzi; ca += cnzi; 259d50806bdSBarry Smith } 260c58ca2e3SHong Zhang 261716bacf3SKris Buschelman ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 262716bacf3SKris Buschelman ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 263d50806bdSBarry Smith ierr = PetscLogFlops(flops);CHKERRQ(ierr); 264d50806bdSBarry Smith PetscFunctionReturn(0); 265d50806bdSBarry Smith } 266bc011b1eSHong Zhang 267e9e4536cSHong Zhang #undef __FUNCT__ 2683c50cad2SHong Zhang #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast" 2693c50cad2SHong Zhang PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast(Mat A,Mat B,PetscReal fill,Mat *C) 27025296bd5SBarry Smith { 27125296bd5SBarry Smith PetscErrorCode ierr; 27225296bd5SBarry Smith Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c; 27325296bd5SBarry Smith PetscInt *ai=a->i,*bi=b->i,*ci,*cj; 2743c50cad2SHong Zhang PetscInt am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N; 27525296bd5SBarry Smith MatScalar *ca; 27625296bd5SBarry Smith PetscReal afill; 2773c50cad2SHong Zhang PetscInt i,j,anzi,brow,bnzj,cnzi,*bj,*aj,nlnk_max,*lnk,ndouble=0; 2783c50cad2SHong Zhang PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 27925296bd5SBarry Smith 28025296bd5SBarry Smith PetscFunctionBegin; 2813c50cad2SHong Zhang /* Get ci and cj - same as MatMatMultSymbolic_SeqAIJ_SeqAIJ except using PetscLLxxx_fast() */ 2823c50cad2SHong Zhang /*-----------------------------------------------------------------------------------------*/ 2833c50cad2SHong Zhang /* Allocate arrays for fill computation and free space for accumulating nonzero column */ 2843c50cad2SHong Zhang ierr = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr); 2853c50cad2SHong Zhang ci[0] = 0; 2863c50cad2SHong Zhang 2873c50cad2SHong Zhang /* create and initialize a linked list */ 2883c50cad2SHong Zhang nlnk_max = a->rmax*b->rmax; 2893c50cad2SHong Zhang if (!nlnk_max || nlnk_max > bn) nlnk_max = bn; /* in case rmax is not defined for A or B */ 2903c50cad2SHong Zhang ierr = PetscLLCondensedCreate_fast(nlnk_max,&lnk);CHKERRQ(ierr); 2913c50cad2SHong Zhang 2923c50cad2SHong Zhang /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */ 2933c50cad2SHong Zhang ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr); 2943c50cad2SHong Zhang current_space = free_space; 2953c50cad2SHong Zhang 2963c50cad2SHong Zhang /* Determine ci and cj */ 2973c50cad2SHong Zhang for (i=0; i<am; i++) { 2983c50cad2SHong Zhang anzi = ai[i+1] - ai[i]; 2993c50cad2SHong Zhang aj = a->j + ai[i]; 3003c50cad2SHong Zhang for (j=0; j<anzi; j++){ 3013c50cad2SHong Zhang brow = aj[j]; 3023c50cad2SHong Zhang bnzj = bi[brow+1] - bi[brow]; 3033c50cad2SHong Zhang bj = b->j + bi[brow]; 3043c50cad2SHong Zhang /* add non-zero cols of B into the sorted linked list lnk */ 3053c50cad2SHong Zhang ierr = PetscLLCondensedAddSorted_fast(bnzj,bj,lnk);CHKERRQ(ierr); 3063c50cad2SHong Zhang } 3073c50cad2SHong Zhang cnzi = lnk[1]; 3083c50cad2SHong Zhang 3093c50cad2SHong Zhang /* If free space is not available, make more free space */ 3103c50cad2SHong Zhang /* Double the amount of total space in the list */ 3113c50cad2SHong Zhang if (current_space->local_remaining<cnzi) { 3123c50cad2SHong Zhang ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 3133c50cad2SHong Zhang ndouble++; 3143c50cad2SHong Zhang } 3153c50cad2SHong Zhang 3163c50cad2SHong Zhang /* Copy data into free space, then initialize lnk */ 3173c50cad2SHong Zhang ierr = PetscLLCondensedClean_fast(cnzi,current_space->array,lnk);CHKERRQ(ierr); 3183c50cad2SHong Zhang current_space->array += cnzi; 3193c50cad2SHong Zhang current_space->local_used += cnzi; 3203c50cad2SHong Zhang current_space->local_remaining -= cnzi; 3213c50cad2SHong Zhang ci[i+1] = ci[i] + cnzi; 3223c50cad2SHong Zhang } 3233c50cad2SHong Zhang 3243c50cad2SHong Zhang /* Column indices are in the list of free space */ 3253c50cad2SHong Zhang /* Allocate space for cj, initialize cj, and */ 3263c50cad2SHong Zhang /* destroy list of free space and other temporary array(s) */ 3273c50cad2SHong Zhang ierr = PetscMalloc((ci[am]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr); 3283c50cad2SHong Zhang ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr); 3293c50cad2SHong Zhang ierr = PetscLLCondensedDestroy_fast(lnk);CHKERRQ(ierr); 33025296bd5SBarry Smith 33125296bd5SBarry Smith /* Allocate space for ca */ 33225296bd5SBarry Smith ierr = PetscMalloc((ci[am]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr); 33325296bd5SBarry Smith ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr); 33425296bd5SBarry Smith 33525296bd5SBarry Smith /* put together the new symbolic matrix */ 33625296bd5SBarry Smith ierr = MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,am,bn,ci,cj,ca,C);CHKERRQ(ierr); 337a2f3521dSMark F. Adams (*C)->rmap->bs = A->rmap->bs; 338a2f3521dSMark F. Adams (*C)->cmap->bs = B->cmap->bs; 33925296bd5SBarry Smith 34025296bd5SBarry Smith /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 34125296bd5SBarry Smith /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */ 34225296bd5SBarry Smith c = (Mat_SeqAIJ *)((*C)->data); 34325296bd5SBarry Smith c->free_a = PETSC_TRUE; 34425296bd5SBarry Smith c->free_ij = PETSC_TRUE; 34525296bd5SBarry Smith c->nonew = 0; 34625296bd5SBarry Smith (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable; /* slower, less memory */ 34725296bd5SBarry Smith 34825296bd5SBarry Smith /* set MatInfo */ 34925296bd5SBarry Smith afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5; 35025296bd5SBarry Smith if (afill < 1.0) afill = 1.0; 35125296bd5SBarry Smith c->maxnz = ci[am]; 35225296bd5SBarry Smith c->nz = ci[am]; 3533c50cad2SHong Zhang (*C)->info.mallocs = ndouble; 35425296bd5SBarry Smith (*C)->info.fill_ratio_given = fill; 35525296bd5SBarry Smith (*C)->info.fill_ratio_needed = afill; 35625296bd5SBarry Smith 35725296bd5SBarry Smith #if defined(PETSC_USE_INFO) 35825296bd5SBarry Smith if (ci[am]) { 3593c50cad2SHong Zhang ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",ndouble,fill,afill);CHKERRQ(ierr); 36025296bd5SBarry Smith ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr); 36125296bd5SBarry Smith } else { 36225296bd5SBarry Smith ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr); 36325296bd5SBarry Smith } 36425296bd5SBarry Smith #endif 36525296bd5SBarry Smith PetscFunctionReturn(0); 36625296bd5SBarry Smith } 36725296bd5SBarry Smith 36825296bd5SBarry Smith 36925296bd5SBarry Smith #undef __FUNCT__ 37025023028SHong Zhang #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable" 37125023028SHong Zhang PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(Mat A,Mat B,PetscReal fill,Mat *C) 372e9e4536cSHong Zhang { 373e9e4536cSHong Zhang PetscErrorCode ierr; 374e9e4536cSHong Zhang Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c; 375bf9555e6SHong Zhang PetscInt *ai=a->i,*bi=b->i,*ci,*cj; 37625c41797SHong Zhang PetscInt am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N; 377e9e4536cSHong Zhang MatScalar *ca; 378e9e4536cSHong Zhang PetscReal afill; 379bd958071SHong Zhang PetscInt i,j,anzi,brow,bnzj,cnzi,*bj,*aj,nlnk_max,*lnk,ndouble=0; 380bd958071SHong Zhang PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 381e9e4536cSHong Zhang 382e9e4536cSHong Zhang PetscFunctionBegin; 383bd958071SHong Zhang /* Get ci and cj - same as MatMatMultSymbolic_SeqAIJ_SeqAIJ except using PetscLLxxx_Scalalbe() */ 384bd958071SHong Zhang /*---------------------------------------------------------------------------------------------*/ 385bd958071SHong Zhang /* Allocate arrays for fill computation and free space for accumulating nonzero column */ 386bd958071SHong Zhang ierr = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr); 387bd958071SHong Zhang ci[0] = 0; 388bd958071SHong Zhang 389bd958071SHong Zhang /* create and initialize a linked list */ 390bd958071SHong Zhang nlnk_max = a->rmax*b->rmax; 391bd958071SHong Zhang if (!nlnk_max || nlnk_max > bn) nlnk_max = bn; /* in case rmax is not defined for A or B */ 392bd958071SHong Zhang ierr = PetscLLCondensedCreate_Scalable(nlnk_max,&lnk);CHKERRQ(ierr); 393bd958071SHong Zhang 394bd958071SHong Zhang /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */ 395bd958071SHong Zhang ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr); 396bd958071SHong Zhang current_space = free_space; 397bd958071SHong Zhang 398bd958071SHong Zhang /* Determine ci and cj */ 399bd958071SHong Zhang for (i=0; i<am; i++) { 400bd958071SHong Zhang anzi = ai[i+1] - ai[i]; 401bd958071SHong Zhang aj = a->j + ai[i]; 402bd958071SHong Zhang for (j=0; j<anzi; j++){ 403bd958071SHong Zhang brow = aj[j]; 404bd958071SHong Zhang bnzj = bi[brow+1] - bi[brow]; 405bd958071SHong Zhang bj = b->j + bi[brow]; 406bd958071SHong Zhang /* add non-zero cols of B into the sorted linked list lnk */ 407bd958071SHong Zhang ierr = PetscLLCondensedAddSorted_Scalable(bnzj,bj,lnk);CHKERRQ(ierr); 408bd958071SHong Zhang } 409bd958071SHong Zhang cnzi = lnk[0]; 410bd958071SHong Zhang 411bd958071SHong Zhang /* If free space is not available, make more free space */ 412bd958071SHong Zhang /* Double the amount of total space in the list */ 413bd958071SHong Zhang if (current_space->local_remaining<cnzi) { 414bd958071SHong Zhang ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 415bd958071SHong Zhang ndouble++; 416bd958071SHong Zhang } 417bd958071SHong Zhang 418bd958071SHong Zhang /* Copy data into free space, then initialize lnk */ 419bd958071SHong Zhang ierr = PetscLLCondensedClean_Scalable(cnzi,current_space->array,lnk);CHKERRQ(ierr); 420bd958071SHong Zhang current_space->array += cnzi; 421bd958071SHong Zhang current_space->local_used += cnzi; 422bd958071SHong Zhang current_space->local_remaining -= cnzi; 423bd958071SHong Zhang ci[i+1] = ci[i] + cnzi; 424bd958071SHong Zhang } 425bd958071SHong Zhang 426bd958071SHong Zhang /* Column indices are in the list of free space */ 427bd958071SHong Zhang /* Allocate space for cj, initialize cj, and */ 428bd958071SHong Zhang /* destroy list of free space and other temporary array(s) */ 429bd958071SHong Zhang ierr = PetscMalloc((ci[am]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr); 430bd958071SHong Zhang ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr); 431bd958071SHong Zhang ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr); 432e9e4536cSHong Zhang 433e9e4536cSHong Zhang /* Allocate space for ca */ 434bd958071SHong Zhang /*-----------------------*/ 435e9e4536cSHong Zhang ierr = PetscMalloc((ci[am]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr); 436e9e4536cSHong Zhang ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr); 437e9e4536cSHong Zhang 438e9e4536cSHong Zhang /* put together the new symbolic matrix */ 439e9e4536cSHong Zhang ierr = MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,am,bn,ci,cj,ca,C);CHKERRQ(ierr); 440a2f3521dSMark F. Adams (*C)->rmap->bs = A->rmap->bs; 441a2f3521dSMark F. Adams (*C)->cmap->bs = B->cmap->bs; 442e9e4536cSHong Zhang 443e9e4536cSHong Zhang /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 444e9e4536cSHong Zhang /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */ 445e9e4536cSHong Zhang c = (Mat_SeqAIJ *)((*C)->data); 446e9e4536cSHong Zhang c->free_a = PETSC_TRUE; 447e9e4536cSHong Zhang c->free_ij = PETSC_TRUE; 448e9e4536cSHong Zhang c->nonew = 0; 44925023028SHong Zhang (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable; /* slower, less memory */ 450e9e4536cSHong Zhang 451e9e4536cSHong Zhang /* set MatInfo */ 452e9e4536cSHong Zhang afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5; 453e9e4536cSHong Zhang if (afill < 1.0) afill = 1.0; 454e9e4536cSHong Zhang c->maxnz = ci[am]; 455e9e4536cSHong Zhang c->nz = ci[am]; 456bd958071SHong Zhang (*C)->info.mallocs = ndouble; 457e9e4536cSHong Zhang (*C)->info.fill_ratio_given = fill; 458e9e4536cSHong Zhang (*C)->info.fill_ratio_needed = afill; 459e9e4536cSHong Zhang 460e9e4536cSHong Zhang #if defined(PETSC_USE_INFO) 461e9e4536cSHong Zhang if (ci[am]) { 462bd958071SHong Zhang ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",ndouble,fill,afill);CHKERRQ(ierr); 463e9e4536cSHong Zhang ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr); 464e9e4536cSHong Zhang } else { 465e9e4536cSHong Zhang ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr); 466e9e4536cSHong Zhang } 467e9e4536cSHong Zhang #endif 468e9e4536cSHong Zhang PetscFunctionReturn(0); 469e9e4536cSHong Zhang } 470e9e4536cSHong Zhang 4710ced3a2bSJed Brown #undef __FUNCT__ 4720ced3a2bSJed Brown #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ_Heap" 4730ced3a2bSJed Brown PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Heap(Mat A,Mat B,PetscReal fill,Mat *C) 4740ced3a2bSJed Brown { 4750ced3a2bSJed Brown PetscErrorCode ierr; 4760ced3a2bSJed Brown Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c; 4770ced3a2bSJed Brown const PetscInt *ai=a->i,*bi=b->i,*aj=a->j,*bj=b->j; 4780ced3a2bSJed Brown PetscInt *ci,*cj,*bb; 4790ced3a2bSJed Brown PetscInt am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N; 4800ced3a2bSJed Brown PetscReal afill; 4810ced3a2bSJed Brown PetscInt i,j,col,ndouble = 0; 4820ced3a2bSJed Brown PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 4830ced3a2bSJed Brown PetscHeap h; 4840ced3a2bSJed Brown 4850ced3a2bSJed Brown PetscFunctionBegin; 4860ced3a2bSJed Brown /* Get ci and cj - same as MatMatMultSymbolic_SeqAIJ_SeqAIJ except using PetscLLxxx_Scalalbe() */ 4870ced3a2bSJed Brown /*---------------------------------------------------------------------------------------------*/ 4880ced3a2bSJed Brown /* Allocate arrays for fill computation and free space for accumulating nonzero column */ 4890ced3a2bSJed Brown ierr = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr); 4900ced3a2bSJed Brown ci[0] = 0; 4910ced3a2bSJed Brown 4920ced3a2bSJed Brown /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */ 4930ced3a2bSJed Brown ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr); 4940ced3a2bSJed Brown current_space = free_space; 4950ced3a2bSJed Brown 4960ced3a2bSJed Brown ierr = PetscHeapCreate(a->rmax,&h);CHKERRQ(ierr); 4970ced3a2bSJed Brown ierr = PetscMalloc(a->rmax*sizeof(PetscInt),&bb);CHKERRQ(ierr); 4980ced3a2bSJed Brown 4990ced3a2bSJed Brown /* Determine ci and cj */ 5000ced3a2bSJed Brown for (i=0; i<am; i++) { 5010ced3a2bSJed 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 */ 5020ced3a2bSJed Brown const PetscInt *acol = aj + ai[i]; /* column indices of nonzero entries in this row */ 5030ced3a2bSJed Brown ci[i+1] = ci[i]; 5040ced3a2bSJed Brown /* Populate the min heap */ 5050ced3a2bSJed Brown for (j=0; j<anzi; j++) { 5060ced3a2bSJed Brown bb[j] = bi[acol[j]]; /* bb points at the start of the row */ 5070ced3a2bSJed Brown if (bb[j] < bi[acol[j]+1]) { /* Add if row is nonempty */ 5080ced3a2bSJed Brown ierr = PetscHeapAdd(h,j,bj[bb[j]++]);CHKERRQ(ierr); 5090ced3a2bSJed Brown } 5100ced3a2bSJed Brown } 5110ced3a2bSJed Brown /* Pick off the min element, adding it to free space */ 5120ced3a2bSJed Brown ierr = PetscHeapPop(h,&j,&col);CHKERRQ(ierr); 5130ced3a2bSJed Brown while (j >= 0) { 5140ced3a2bSJed Brown if (current_space->local_remaining < 1) { /* double the size, but don't exceed 16 MiB */ 5150ced3a2bSJed Brown ierr = PetscFreeSpaceGet(PetscMin(2*current_space->total_array_size,16 << 20),¤t_space);CHKERRQ(ierr); 5160ced3a2bSJed Brown ndouble++; 5170ced3a2bSJed Brown } 5180ced3a2bSJed Brown *(current_space->array++) = col; 5190ced3a2bSJed Brown current_space->local_used++; 5200ced3a2bSJed Brown current_space->local_remaining--; 5210ced3a2bSJed Brown ci[i+1]++; 5220ced3a2bSJed Brown 5230ced3a2bSJed Brown /* stash if anything else remains in this row of B */ 5240ced3a2bSJed Brown if (bb[j] < bi[acol[j]+1]) {ierr = PetscHeapStash(h,j,bj[bb[j]++]);CHKERRQ(ierr);} 5250ced3a2bSJed Brown while (1) { /* pop and stash any other rows of B that also had an entry in this column */ 5260ced3a2bSJed Brown PetscInt j2,col2; 5270ced3a2bSJed Brown ierr = PetscHeapPeek(h,&j2,&col2);CHKERRQ(ierr); 5280ced3a2bSJed Brown if (col2 != col) break; 5290ced3a2bSJed Brown ierr = PetscHeapPop(h,&j2,&col2);CHKERRQ(ierr); 5300ced3a2bSJed Brown if (bb[j2] < bi[acol[j2]+1]) {ierr = PetscHeapStash(h,j2,bj[bb[j2]++]);CHKERRQ(ierr);} 5310ced3a2bSJed Brown } 5320ced3a2bSJed Brown /* Put any stashed elements back into the min heap */ 5330ced3a2bSJed Brown ierr = PetscHeapUnstash(h);CHKERRQ(ierr); 5340ced3a2bSJed Brown ierr = PetscHeapPop(h,&j,&col);CHKERRQ(ierr); 5350ced3a2bSJed Brown } 5360ced3a2bSJed Brown } 5370ced3a2bSJed Brown ierr = PetscFree(bb);CHKERRQ(ierr); 5380ced3a2bSJed Brown ierr = PetscHeapDestroy(&h);CHKERRQ(ierr); 5390ced3a2bSJed Brown 5400ced3a2bSJed Brown /* Column indices are in the list of free space */ 5410ced3a2bSJed Brown /* Allocate space for cj, initialize cj, and */ 5420ced3a2bSJed Brown /* destroy list of free space and other temporary array(s) */ 5430ced3a2bSJed Brown ierr = PetscMalloc(ci[am]*sizeof(PetscInt),&cj);CHKERRQ(ierr); 5440ced3a2bSJed Brown ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr); 5450ced3a2bSJed Brown 5460ced3a2bSJed Brown /* put together the new symbolic matrix */ 54789d95c1aSJed Brown ierr = MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,am,bn,ci,cj,PETSC_NULL,C);CHKERRQ(ierr); 5480ced3a2bSJed Brown (*C)->rmap->bs = A->rmap->bs; 5490ced3a2bSJed Brown (*C)->cmap->bs = B->cmap->bs; 5500ced3a2bSJed Brown 5510ced3a2bSJed Brown /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 5520ced3a2bSJed Brown /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */ 5530ced3a2bSJed Brown c = (Mat_SeqAIJ *)((*C)->data); 5540ced3a2bSJed Brown c->free_a = PETSC_TRUE; 5550ced3a2bSJed Brown c->free_ij = PETSC_TRUE; 5560ced3a2bSJed Brown c->nonew = 0; 55789d95c1aSJed Brown (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ; 5580ced3a2bSJed Brown 5590ced3a2bSJed Brown /* set MatInfo */ 5600ced3a2bSJed Brown afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5; 5610ced3a2bSJed Brown if (afill < 1.0) afill = 1.0; 5620ced3a2bSJed Brown c->maxnz = ci[am]; 5630ced3a2bSJed Brown c->nz = ci[am]; 5640ced3a2bSJed Brown (*C)->info.mallocs = ndouble; 5650ced3a2bSJed Brown (*C)->info.fill_ratio_given = fill; 5660ced3a2bSJed Brown (*C)->info.fill_ratio_needed = afill; 5670ced3a2bSJed Brown 5680ced3a2bSJed Brown #if defined(PETSC_USE_INFO) 5690ced3a2bSJed Brown if (ci[am]) { 5700ced3a2bSJed Brown ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",ndouble,fill,afill);CHKERRQ(ierr); 5710ced3a2bSJed Brown ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr); 5720ced3a2bSJed Brown } else { 5730ced3a2bSJed Brown ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr); 5740ced3a2bSJed Brown } 5750ced3a2bSJed Brown #endif 5760ced3a2bSJed Brown PetscFunctionReturn(0); 5770ced3a2bSJed Brown } 578e9e4536cSHong Zhang 5798a07c6f1SJed Brown #undef __FUNCT__ 5808a07c6f1SJed Brown #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ_BTHeap" 5818a07c6f1SJed Brown PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_BTHeap(Mat A,Mat B,PetscReal fill,Mat *C) 5828a07c6f1SJed Brown { 5838a07c6f1SJed Brown PetscErrorCode ierr; 5848a07c6f1SJed Brown Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c; 5858a07c6f1SJed Brown const PetscInt *ai=a->i,*bi=b->i,*aj=a->j,*bj=b->j; 5868a07c6f1SJed Brown PetscInt *ci,*cj,*bb; 5878a07c6f1SJed Brown PetscInt am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N; 5888a07c6f1SJed Brown PetscReal afill; 5898a07c6f1SJed Brown PetscInt i,j,col,ndouble = 0; 5908a07c6f1SJed Brown PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 5918a07c6f1SJed Brown PetscHeap h; 5928a07c6f1SJed Brown PetscBT bt; 5938a07c6f1SJed Brown 5948a07c6f1SJed Brown PetscFunctionBegin; 5958a07c6f1SJed Brown /* Get ci and cj - same as MatMatMultSymbolic_SeqAIJ_SeqAIJ except using PetscLLxxx_Scalalbe() */ 5968a07c6f1SJed Brown /*---------------------------------------------------------------------------------------------*/ 5978a07c6f1SJed Brown /* Allocate arrays for fill computation and free space for accumulating nonzero column */ 5988a07c6f1SJed Brown ierr = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr); 5998a07c6f1SJed Brown ci[0] = 0; 6008a07c6f1SJed Brown 6018a07c6f1SJed Brown /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */ 6028a07c6f1SJed Brown ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr); 6038a07c6f1SJed Brown current_space = free_space; 6048a07c6f1SJed Brown 6058a07c6f1SJed Brown ierr = PetscHeapCreate(a->rmax,&h);CHKERRQ(ierr); 6068a07c6f1SJed Brown ierr = PetscMalloc(a->rmax*sizeof(PetscInt),&bb);CHKERRQ(ierr); 6078a07c6f1SJed Brown ierr = PetscBTCreate(bn,&bt);CHKERRQ(ierr); 6088a07c6f1SJed Brown 6098a07c6f1SJed Brown /* Determine ci and cj */ 6108a07c6f1SJed Brown for (i=0; i<am; i++) { 6118a07c6f1SJed 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 */ 6128a07c6f1SJed Brown const PetscInt *acol = aj + ai[i]; /* column indices of nonzero entries in this row */ 6138a07c6f1SJed Brown const PetscInt *fptr = current_space->array; /* Save beginning of the row so we can clear the BT later */ 6148a07c6f1SJed Brown ci[i+1] = ci[i]; 6158a07c6f1SJed Brown /* Populate the min heap */ 6168a07c6f1SJed Brown for (j=0; j<anzi; j++) { 6178a07c6f1SJed Brown PetscInt brow = acol[j]; 6188a07c6f1SJed Brown for (bb[j] = bi[brow]; bb[j] < bi[brow+1]; bb[j]++) { 6198a07c6f1SJed Brown PetscInt bcol = bj[bb[j]]; 6208a07c6f1SJed Brown if (!PetscBTLookupSet(bt,bcol)) { /* new entry */ 6218a07c6f1SJed Brown ierr = PetscHeapAdd(h,j,bcol);CHKERRQ(ierr); 6228a07c6f1SJed Brown bb[j]++; 6238a07c6f1SJed Brown break; 6248a07c6f1SJed Brown } 6258a07c6f1SJed Brown } 6268a07c6f1SJed Brown } 6278a07c6f1SJed Brown /* Pick off the min element, adding it to free space */ 6288a07c6f1SJed Brown ierr = PetscHeapPop(h,&j,&col);CHKERRQ(ierr); 6298a07c6f1SJed Brown while (j >= 0) { 6308a07c6f1SJed Brown if (current_space->local_remaining < 1) { /* double the size, but don't exceed 16 MiB */ 6318a07c6f1SJed Brown fptr = PETSC_NULL; /* need PetscBTMemzero */ 6328a07c6f1SJed Brown ierr = PetscFreeSpaceGet(PetscMin(2*current_space->total_array_size,16 << 20),¤t_space);CHKERRQ(ierr); 6338a07c6f1SJed Brown ndouble++; 6348a07c6f1SJed Brown } 6358a07c6f1SJed Brown *(current_space->array++) = col; 6368a07c6f1SJed Brown current_space->local_used++; 6378a07c6f1SJed Brown current_space->local_remaining--; 6388a07c6f1SJed Brown ci[i+1]++; 6398a07c6f1SJed Brown 6408a07c6f1SJed Brown /* stash if anything else remains in this row of B */ 6418a07c6f1SJed Brown for ( ; bb[j] < bi[acol[j]+1]; bb[j]++) { 6428a07c6f1SJed Brown PetscInt bcol = bj[bb[j]]; 6438a07c6f1SJed Brown if (!PetscBTLookupSet(bt,bcol)) { /* new entry */ 6448a07c6f1SJed Brown ierr = PetscHeapAdd(h,j,bcol);CHKERRQ(ierr); 6458a07c6f1SJed Brown bb[j]++; 6468a07c6f1SJed Brown break; 6478a07c6f1SJed Brown } 6488a07c6f1SJed Brown } 6498a07c6f1SJed Brown ierr = PetscHeapPop(h,&j,&col);CHKERRQ(ierr); 6508a07c6f1SJed Brown } 6518a07c6f1SJed Brown if (fptr) { /* Clear the bits for this row */ 6528a07c6f1SJed Brown for ( ; fptr<current_space->array; fptr++) {ierr = PetscBTClear(bt,*fptr);CHKERRQ(ierr);} 6538a07c6f1SJed Brown } else { /* We reallocated so we don't remember (easily) how to clear only the bits we changed */ 6548a07c6f1SJed Brown ierr = PetscBTMemzero(bn,bt);CHKERRQ(ierr); 6558a07c6f1SJed Brown } 6568a07c6f1SJed Brown } 6578a07c6f1SJed Brown ierr = PetscFree(bb);CHKERRQ(ierr); 6588a07c6f1SJed Brown ierr = PetscHeapDestroy(&h);CHKERRQ(ierr); 6598a07c6f1SJed Brown ierr = PetscBTDestroy(&bt);CHKERRQ(ierr); 6608a07c6f1SJed Brown 6618a07c6f1SJed Brown /* Column indices are in the list of free space */ 6628a07c6f1SJed Brown /* Allocate space for cj, initialize cj, and */ 6638a07c6f1SJed Brown /* destroy list of free space and other temporary array(s) */ 6648a07c6f1SJed Brown ierr = PetscMalloc(ci[am]*sizeof(PetscInt),&cj);CHKERRQ(ierr); 6658a07c6f1SJed Brown ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr); 6668a07c6f1SJed Brown 6678a07c6f1SJed Brown /* put together the new symbolic matrix */ 66889d95c1aSJed Brown ierr = MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,am,bn,ci,cj,PETSC_NULL,C);CHKERRQ(ierr); 6698a07c6f1SJed Brown (*C)->rmap->bs = A->rmap->bs; 6708a07c6f1SJed Brown (*C)->cmap->bs = B->cmap->bs; 6718a07c6f1SJed Brown 6728a07c6f1SJed Brown /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 6738a07c6f1SJed Brown /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */ 6748a07c6f1SJed Brown c = (Mat_SeqAIJ *)((*C)->data); 6758a07c6f1SJed Brown c->free_a = PETSC_TRUE; 6768a07c6f1SJed Brown c->free_ij = PETSC_TRUE; 6778a07c6f1SJed Brown c->nonew = 0; 67889d95c1aSJed Brown (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ; 6798a07c6f1SJed Brown 6808a07c6f1SJed Brown /* set MatInfo */ 6818a07c6f1SJed Brown afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5; 6828a07c6f1SJed Brown if (afill < 1.0) afill = 1.0; 6838a07c6f1SJed Brown c->maxnz = ci[am]; 6848a07c6f1SJed Brown c->nz = ci[am]; 6858a07c6f1SJed Brown (*C)->info.mallocs = ndouble; 6868a07c6f1SJed Brown (*C)->info.fill_ratio_given = fill; 6878a07c6f1SJed Brown (*C)->info.fill_ratio_needed = afill; 6888a07c6f1SJed Brown 6898a07c6f1SJed Brown #if defined(PETSC_USE_INFO) 6908a07c6f1SJed Brown if (ci[am]) { 6918a07c6f1SJed Brown ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",ndouble,fill,afill);CHKERRQ(ierr); 6928a07c6f1SJed Brown ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr); 6938a07c6f1SJed Brown } else { 6948a07c6f1SJed Brown ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr); 6958a07c6f1SJed Brown } 6968a07c6f1SJed Brown #endif 6978a07c6f1SJed Brown PetscFunctionReturn(0); 6988a07c6f1SJed Brown } 6998a07c6f1SJed Brown 700d2853540SHong Zhang /* This routine is not used. Should be removed! */ 701bc011b1eSHong Zhang #undef __FUNCT__ 7026fc122caSHong Zhang #define __FUNCT__ "MatMatTransposeMult_SeqAIJ_SeqAIJ" 7036fc122caSHong Zhang PetscErrorCode MatMatTransposeMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 7045df89d91SHong Zhang { 705bc011b1eSHong Zhang PetscErrorCode ierr; 706bc011b1eSHong Zhang 707bc011b1eSHong Zhang PetscFunctionBegin; 708bc011b1eSHong Zhang if (scall == MAT_INITIAL_MATRIX){ 7096fc122caSHong Zhang ierr = MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr); 710bc011b1eSHong Zhang } 7116fc122caSHong Zhang ierr = MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr); 712bc011b1eSHong Zhang PetscFunctionReturn(0); 713bc011b1eSHong Zhang } 714bc011b1eSHong Zhang 715bc011b1eSHong Zhang #undef __FUNCT__ 716e286af10SHong Zhang #define __FUNCT__ "PetscContainerDestroy_Mat_MatMatTransMult" 717e286af10SHong Zhang PetscErrorCode PetscContainerDestroy_Mat_MatMatTransMult(void *ptr) 7182128a86cSHong Zhang { 7192128a86cSHong Zhang PetscErrorCode ierr; 720e286af10SHong Zhang Mat_MatMatTransMult *multtrans=(Mat_MatMatTransMult*)ptr; 7212128a86cSHong Zhang 7222128a86cSHong Zhang PetscFunctionBegin; 723b9af6bddSHong Zhang ierr = MatTransposeColoringDestroy(&multtrans->matcoloring);CHKERRQ(ierr); 7242128a86cSHong Zhang ierr = MatDestroy(&multtrans->Bt_den);CHKERRQ(ierr); 7252128a86cSHong Zhang ierr = MatDestroy(&multtrans->ABt_den);CHKERRQ(ierr); 7262128a86cSHong Zhang ierr = PetscFree(multtrans);CHKERRQ(ierr); 7272128a86cSHong Zhang PetscFunctionReturn(0); 7282128a86cSHong Zhang } 7292128a86cSHong Zhang 7302128a86cSHong Zhang #undef __FUNCT__ 7312128a86cSHong Zhang #define __FUNCT__ "MatDestroy_SeqAIJ_MatMatMultTrans" 7322128a86cSHong Zhang PetscErrorCode MatDestroy_SeqAIJ_MatMatMultTrans(Mat A) 7332128a86cSHong Zhang { 7342128a86cSHong Zhang PetscErrorCode ierr; 7352128a86cSHong Zhang PetscContainer container; 736e286af10SHong Zhang Mat_MatMatTransMult *multtrans=PETSC_NULL; 7372128a86cSHong Zhang 7382128a86cSHong Zhang PetscFunctionBegin; 739e286af10SHong Zhang ierr = PetscObjectQuery((PetscObject)A,"Mat_MatMatTransMult",(PetscObject *)&container);CHKERRQ(ierr); 7402128a86cSHong Zhang if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit"); 7412128a86cSHong Zhang ierr = PetscContainerGetPointer(container,(void **)&multtrans);CHKERRQ(ierr); 7422128a86cSHong Zhang A->ops->destroy = multtrans->destroy; 7432128a86cSHong Zhang if (A->ops->destroy) { 7442128a86cSHong Zhang ierr = (*A->ops->destroy)(A);CHKERRQ(ierr); 7452128a86cSHong Zhang } 746e286af10SHong Zhang ierr = PetscObjectCompose((PetscObject)A,"Mat_MatMatTransMult",0);CHKERRQ(ierr); 7472128a86cSHong Zhang PetscFunctionReturn(0); 7482128a86cSHong Zhang } 7492128a86cSHong Zhang 7502128a86cSHong Zhang #undef __FUNCT__ 7516fc122caSHong Zhang #define __FUNCT__ "MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ" 7526fc122caSHong Zhang PetscErrorCode MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C) 753bc011b1eSHong Zhang { 7545df89d91SHong Zhang PetscErrorCode ierr; 75581d82fe4SHong Zhang Mat Bt; 75681d82fe4SHong Zhang PetscInt *bti,*btj; 757e286af10SHong Zhang Mat_MatMatTransMult *multtrans; 7582128a86cSHong Zhang PetscContainer container; 759d2853540SHong Zhang PetscLogDouble t0,tf,etime2=0.0; 760d2853540SHong Zhang 76181d82fe4SHong Zhang PetscFunctionBegin; 762d2853540SHong Zhang ierr = PetscGetTime(&t0);CHKERRQ(ierr); 76381d82fe4SHong Zhang /* create symbolic Bt */ 76481d82fe4SHong Zhang ierr = MatGetSymbolicTranspose_SeqAIJ(B,&bti,&btj);CHKERRQ(ierr); 76581d82fe4SHong Zhang ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,B->cmap->n,B->rmap->n,bti,btj,PETSC_NULL,&Bt);CHKERRQ(ierr); 766c3e5699bSSatish Balay Bt->rmap->bs = A->cmap->bs; 767c3e5699bSSatish Balay Bt->cmap->bs = B->cmap->bs; 76881d82fe4SHong Zhang 76981d82fe4SHong Zhang /* get symbolic C=A*Bt */ 77081d82fe4SHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,Bt,fill,C);CHKERRQ(ierr); 77181d82fe4SHong Zhang 7722128a86cSHong Zhang /* create a supporting struct for reuse intermidiate dense matrices with matcoloring */ 773e286af10SHong Zhang ierr = PetscNew(Mat_MatMatTransMult,&multtrans);CHKERRQ(ierr); 7742128a86cSHong Zhang 7752128a86cSHong Zhang /* attach the supporting struct to C */ 7762128a86cSHong Zhang ierr = PetscContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr); 7772128a86cSHong Zhang ierr = PetscContainerSetPointer(container,multtrans);CHKERRQ(ierr); 778e286af10SHong Zhang ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_Mat_MatMatTransMult);CHKERRQ(ierr); 779e286af10SHong Zhang ierr = PetscObjectCompose((PetscObject)(*C),"Mat_MatMatTransMult",(PetscObject)container);CHKERRQ(ierr); 7802128a86cSHong Zhang ierr = PetscContainerDestroy(&container);CHKERRQ(ierr); 7812128a86cSHong Zhang 7822128a86cSHong Zhang multtrans->usecoloring = PETSC_FALSE; 7832128a86cSHong Zhang multtrans->destroy = (*C)->ops->destroy; 7842128a86cSHong Zhang (*C)->ops->destroy = MatDestroy_SeqAIJ_MatMatMultTrans; 7852128a86cSHong Zhang 786d2853540SHong Zhang ierr = PetscGetTime(&tf);CHKERRQ(ierr); 787d2853540SHong Zhang etime2 += tf - t0; 788d2853540SHong Zhang 789f400d4cbSHong Zhang ierr = PetscOptionsGetBool(PETSC_NULL,"-matmattransmult_color",&multtrans->usecoloring,PETSC_NULL);CHKERRQ(ierr); 7902128a86cSHong Zhang if (multtrans->usecoloring){ 791b9af6bddSHong Zhang /* Create MatTransposeColoring from symbolic C=A*B^T */ 792b9af6bddSHong Zhang MatTransposeColoring matcoloring; 7932128a86cSHong Zhang ISColoring iscoloring; 7942128a86cSHong Zhang Mat Bt_dense,C_dense; 795d2853540SHong Zhang PetscLogDouble etime0=0.0,etime01=0.0,etime1=0.0; 796e8354b3eSHong Zhang 797e8354b3eSHong Zhang ierr = PetscGetTime(&t0);CHKERRQ(ierr); 798502de53fSHong Zhang ierr = MatGetColoring(*C,MATCOLORINGLF,&iscoloring);CHKERRQ(ierr); 799d2853540SHong Zhang ierr = PetscGetTime(&tf);CHKERRQ(ierr); 800d2853540SHong Zhang etime0 += tf - t0; 801d2853540SHong Zhang 802d2853540SHong Zhang ierr = PetscGetTime(&t0);CHKERRQ(ierr); 803b9af6bddSHong Zhang ierr = MatTransposeColoringCreate(*C,iscoloring,&matcoloring);CHKERRQ(ierr); 8042128a86cSHong Zhang multtrans->matcoloring = matcoloring; 8052128a86cSHong Zhang ierr = ISColoringDestroy(&iscoloring);CHKERRQ(ierr); 806e8354b3eSHong Zhang ierr = PetscGetTime(&tf);CHKERRQ(ierr); 807d2853540SHong Zhang etime01 += tf - t0; 8082128a86cSHong Zhang 809e8354b3eSHong Zhang ierr = PetscGetTime(&t0);CHKERRQ(ierr); 8102128a86cSHong Zhang /* Create Bt_dense and C_dense = A*Bt_dense */ 8112128a86cSHong Zhang ierr = MatCreate(PETSC_COMM_SELF,&Bt_dense);CHKERRQ(ierr); 8122128a86cSHong Zhang ierr = MatSetSizes(Bt_dense,A->cmap->n,matcoloring->ncolors,A->cmap->n,matcoloring->ncolors);CHKERRQ(ierr); 8132128a86cSHong Zhang ierr = MatSetType(Bt_dense,MATSEQDENSE);CHKERRQ(ierr); 8142128a86cSHong Zhang ierr = MatSeqDenseSetPreallocation(Bt_dense,PETSC_NULL);CHKERRQ(ierr); 8152128a86cSHong Zhang Bt_dense->assembled = PETSC_TRUE; 8162128a86cSHong Zhang multtrans->Bt_den = Bt_dense; 8172128a86cSHong Zhang 8182128a86cSHong Zhang ierr = MatCreate(PETSC_COMM_SELF,&C_dense);CHKERRQ(ierr); 8192128a86cSHong Zhang ierr = MatSetSizes(C_dense,A->rmap->n,matcoloring->ncolors,A->rmap->n,matcoloring->ncolors);CHKERRQ(ierr); 8202128a86cSHong Zhang ierr = MatSetType(C_dense,MATSEQDENSE);CHKERRQ(ierr); 8212128a86cSHong Zhang ierr = MatSeqDenseSetPreallocation(C_dense,PETSC_NULL);CHKERRQ(ierr); 8222128a86cSHong Zhang Bt_dense->assembled = PETSC_TRUE; 8232128a86cSHong Zhang multtrans->ABt_den = C_dense; 824e8354b3eSHong Zhang ierr = PetscGetTime(&tf);CHKERRQ(ierr); 825e8354b3eSHong Zhang etime1 += tf - t0; 826f94ccd6cSHong Zhang 827f94ccd6cSHong Zhang #if defined(PETSC_USE_INFO) 8281ce71dffSSatish Balay { 829f94ccd6cSHong Zhang Mat_SeqAIJ *c=(Mat_SeqAIJ*)(*C)->data; 830f94ccd6cSHong 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)); 831f94ccd6cSHong Zhang ierr = PetscInfo5(*C,"Sym = GetColor %g + ColorCreate %g + MatDenseCreate %g + non-colorSym %g = %g\n",etime0,etime01,etime1,etime2,etime0+etime01+etime1+etime2); 8321ce71dffSSatish Balay } 833f94ccd6cSHong Zhang #endif 8342128a86cSHong Zhang } 835e99be685SHong Zhang /* clean up */ 836e99be685SHong Zhang ierr = MatDestroy(&Bt);CHKERRQ(ierr); 837e99be685SHong Zhang ierr = MatRestoreSymbolicTranspose_SeqAIJ(B,&bti,&btj);CHKERRQ(ierr); 8382128a86cSHong Zhang 839f94ccd6cSHong Zhang 840f94ccd6cSHong Zhang 84181d82fe4SHong Zhang #if defined(INEFFICIENT_ALGORITHM) 84281d82fe4SHong Zhang /* The algorithm below computes am*bm sparse inner-product - inefficient! It will be deleted later. */ 8431fa1209cSHong Zhang PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 8441fa1209cSHong Zhang Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c; 8451fa1209cSHong Zhang PetscInt *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*ci,*cj,*acol,*bcol; 8461fa1209cSHong Zhang PetscInt am=A->rmap->N,bm=B->rmap->N; 8471fa1209cSHong Zhang PetscInt i,j,anzi,bnzj,cnzi,nlnk,*lnk,nspacedouble=0,ka,kb,index[1]; 8481fa1209cSHong Zhang MatScalar *ca; 8491fa1209cSHong Zhang PetscBT lnkbt; 8501fa1209cSHong Zhang PetscReal afill; 8515df89d91SHong Zhang 8521fa1209cSHong Zhang /* Allocate row pointer array ci */ 8531fa1209cSHong Zhang ierr = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr); 8541fa1209cSHong Zhang ci[0] = 0; 8551fa1209cSHong Zhang 8561fa1209cSHong Zhang /* Create and initialize a linked list for C columns */ 8571fa1209cSHong Zhang nlnk = bm+1; 8581fa1209cSHong Zhang ierr = PetscLLCreate(bm,bm,nlnk,lnk,lnkbt);CHKERRQ(ierr); 8591fa1209cSHong Zhang 8601fa1209cSHong Zhang /* Initial FreeSpace with size fill*(nnz(A)+nnz(B)) */ 8611fa1209cSHong Zhang ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr); 8621fa1209cSHong Zhang current_space = free_space; 8631fa1209cSHong Zhang 8641fa1209cSHong Zhang /* Determine symbolic info for each row of the product A*B^T: */ 8651fa1209cSHong Zhang for (i=0; i<am; i++) { 8661fa1209cSHong Zhang anzi = ai[i+1] - ai[i]; 8671fa1209cSHong Zhang cnzi = 0; 8681fa1209cSHong Zhang acol = aj + ai[i]; 8691fa1209cSHong Zhang for (j=0; j<bm; j++){ 8701fa1209cSHong Zhang bnzj = bi[j+1] - bi[j]; 8711fa1209cSHong Zhang bcol= bj + bi[j]; 87281d82fe4SHong Zhang /* sparse inner-product c(i,j)=A[i,:]*B[j,:]^T */ 8731fa1209cSHong Zhang ka = 0; kb = 0; 8741fa1209cSHong Zhang while (ka < anzi && kb < bnzj){ 87581d82fe4SHong Zhang while (acol[ka] < bcol[kb] && ka < anzi) ka++; 87681d82fe4SHong Zhang if (ka == anzi) break; 87781d82fe4SHong Zhang while (acol[ka] > bcol[kb] && kb < bnzj) kb++; 87881d82fe4SHong Zhang if (kb == bnzj) break; 8791fa1209cSHong Zhang if (acol[ka] == bcol[kb]){ /* add nonzero c(i,j) to lnk */ 8801fa1209cSHong Zhang index[0] = j; 8811fa1209cSHong Zhang ierr = PetscLLAdd(1,index,bm,nlnk,lnk,lnkbt);CHKERRQ(ierr); 8821fa1209cSHong Zhang cnzi++; 8831fa1209cSHong Zhang break; 8841fa1209cSHong Zhang } 8851fa1209cSHong Zhang } 8861fa1209cSHong Zhang } 8871fa1209cSHong Zhang 8881fa1209cSHong Zhang /* If free space is not available, make more free space */ 8891fa1209cSHong Zhang /* Double the amount of total space in the list */ 8901fa1209cSHong Zhang if (current_space->local_remaining<cnzi) { 8911fa1209cSHong Zhang ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 8921fa1209cSHong Zhang nspacedouble++; 8931fa1209cSHong Zhang } 8941fa1209cSHong Zhang 8951fa1209cSHong Zhang /* Copy data into free space, then initialize lnk */ 8961fa1209cSHong Zhang ierr = PetscLLClean(bm,bm,cnzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 8971fa1209cSHong Zhang current_space->array += cnzi; 8981fa1209cSHong Zhang current_space->local_used += cnzi; 8991fa1209cSHong Zhang current_space->local_remaining -= cnzi; 9001fa1209cSHong Zhang 9011fa1209cSHong Zhang ci[i+1] = ci[i] + cnzi; 9021fa1209cSHong Zhang } 9031fa1209cSHong Zhang 9041fa1209cSHong Zhang 9051fa1209cSHong Zhang /* Column indices are in the list of free space. 9061fa1209cSHong Zhang Allocate array cj, copy column indices to cj, and destroy list of free space */ 9071fa1209cSHong Zhang ierr = PetscMalloc((ci[am]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr); 9081fa1209cSHong Zhang ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr); 9091fa1209cSHong Zhang ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 9101fa1209cSHong Zhang 9111fa1209cSHong Zhang /* Allocate space for ca */ 9121fa1209cSHong Zhang ierr = PetscMalloc((ci[am]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr); 9131fa1209cSHong Zhang ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr); 9141fa1209cSHong Zhang 9151fa1209cSHong Zhang /* put together the new symbolic matrix */ 9161fa1209cSHong Zhang ierr = MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,am,bm,ci,cj,ca,C);CHKERRQ(ierr); 917a2f3521dSMark F. Adams (*C)->rmap->bs = A->cmap->bs; 918a2f3521dSMark F. Adams (*C)->cmap->bs = B->cmap->bs; 9191fa1209cSHong Zhang 9201fa1209cSHong Zhang /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 9211fa1209cSHong Zhang /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */ 9221fa1209cSHong Zhang c = (Mat_SeqAIJ *)((*C)->data); 9231fa1209cSHong Zhang c->free_a = PETSC_TRUE; 9241fa1209cSHong Zhang c->free_ij = PETSC_TRUE; 9251fa1209cSHong Zhang c->nonew = 0; 9261fa1209cSHong Zhang 9271fa1209cSHong Zhang /* set MatInfo */ 9281fa1209cSHong Zhang afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5; 9291fa1209cSHong Zhang if (afill < 1.0) afill = 1.0; 9301fa1209cSHong Zhang c->maxnz = ci[am]; 9311fa1209cSHong Zhang c->nz = ci[am]; 9321fa1209cSHong Zhang (*C)->info.mallocs = nspacedouble; 9331fa1209cSHong Zhang (*C)->info.fill_ratio_given = fill; 9341fa1209cSHong Zhang (*C)->info.fill_ratio_needed = afill; 9351fa1209cSHong Zhang 9361fa1209cSHong Zhang #if defined(PETSC_USE_INFO) 9371fa1209cSHong Zhang if (ci[am]) { 9381fa1209cSHong Zhang ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr); 9396fc122caSHong Zhang ierr = PetscInfo1((*C),"Use MatMatTransposeMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr); 9401fa1209cSHong Zhang } else { 9411fa1209cSHong Zhang ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr); 9421fa1209cSHong Zhang } 9431fa1209cSHong Zhang #endif 94481d82fe4SHong Zhang #endif 9455df89d91SHong Zhang PetscFunctionReturn(0); 9465df89d91SHong Zhang } 9475df89d91SHong Zhang 9486973769bSHong Zhang /* #define USE_ARRAY - for sparse dot product. Slower than !USE_ARRAY */ 9495df89d91SHong Zhang #undef __FUNCT__ 9506fc122caSHong Zhang #define __FUNCT__ "MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ" 9516fc122caSHong Zhang PetscErrorCode MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C) 9525df89d91SHong Zhang { 9535df89d91SHong Zhang PetscErrorCode ierr; 9545df89d91SHong Zhang Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data; 955e2cac8caSJed Brown PetscInt *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,anzi,bnzj,nexta,nextb,*acol,*bcol,brow; 95681d82fe4SHong Zhang PetscInt cm=C->rmap->n,*ci=c->i,*cj=c->j,i,j,cnzi,*ccol; 9575df89d91SHong Zhang PetscLogDouble flops=0.0; 958aa1aec99SHong Zhang MatScalar *aa=a->a,*aval,*ba=b->a,*bval,*ca,*cval; 959e286af10SHong Zhang Mat_MatMatTransMult *multtrans; 9602128a86cSHong Zhang PetscContainer container; 9616973769bSHong Zhang #if defined(USE_ARRAY) 9626973769bSHong Zhang MatScalar *spdot; 9636973769bSHong Zhang #endif 9645df89d91SHong Zhang 9655df89d91SHong Zhang PetscFunctionBegin; 966e286af10SHong Zhang ierr = PetscObjectQuery((PetscObject)C,"Mat_MatMatTransMult",(PetscObject *)&container);CHKERRQ(ierr); 9672128a86cSHong Zhang if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit"); 9682128a86cSHong Zhang ierr = PetscContainerGetPointer(container,(void **)&multtrans);CHKERRQ(ierr); 9692128a86cSHong Zhang if (multtrans->usecoloring){ 970b9af6bddSHong Zhang MatTransposeColoring matcoloring = multtrans->matcoloring; 971c8db22f6SHong Zhang Mat Bt_dense; 972c8db22f6SHong Zhang PetscInt m,n; 973b2d2b04fSHong Zhang PetscLogDouble t0,tf,etime0=0.0,etime1=0.0,etime2=0.0; 974a0b95ffbSSatish Balay Mat C_dense = multtrans->ABt_den; 975a0b95ffbSSatish Balay 9762128a86cSHong Zhang Bt_dense = multtrans->Bt_den; 977c8db22f6SHong Zhang ierr = MatGetLocalSize(Bt_dense,&m,&n);CHKERRQ(ierr); 978c8db22f6SHong Zhang 979b9af6bddSHong Zhang /* Get Bt_dense by Apply MatTransposeColoring to B */ 9802128a86cSHong Zhang ierr = PetscGetTime(&t0);CHKERRQ(ierr); 981b9af6bddSHong Zhang ierr = MatTransColoringApplySpToDen(matcoloring,B,Bt_dense);CHKERRQ(ierr); 9822128a86cSHong Zhang ierr = PetscGetTime(&tf);CHKERRQ(ierr); 983b2d2b04fSHong Zhang etime0 += tf - t0; 984c8db22f6SHong Zhang 985c8db22f6SHong Zhang /* C_dense = A*Bt_dense */ 9862128a86cSHong Zhang ierr = PetscGetTime(&t0);CHKERRQ(ierr); 9872128a86cSHong Zhang ierr = MatMatMultNumeric_SeqAIJ_SeqDense(A,Bt_dense,C_dense);CHKERRQ(ierr); 9882128a86cSHong Zhang ierr = PetscGetTime(&tf);CHKERRQ(ierr); 9892128a86cSHong Zhang etime2 += tf - t0; 990c8db22f6SHong Zhang 9912128a86cSHong Zhang /* Recover C from C_dense */ 9922128a86cSHong Zhang ierr = PetscGetTime(&t0);CHKERRQ(ierr); 993b9af6bddSHong Zhang ierr = MatTransColoringApplyDenToSp(matcoloring,C_dense,C);CHKERRQ(ierr); 9942128a86cSHong Zhang ierr = PetscGetTime(&tf);CHKERRQ(ierr); 9952128a86cSHong Zhang etime1 += tf - t0; 996f94ccd6cSHong Zhang #if defined(PETSC_USE_INFO) 997f94ccd6cSHong Zhang ierr = PetscInfo4(C,"Num = ColoringApply: %g %g + Mult_sp_dense %g = %g\n",etime0,etime1,etime2,etime0+etime1+etime2); 998f94ccd6cSHong Zhang #endif 999c8db22f6SHong Zhang PetscFunctionReturn(0); 1000c8db22f6SHong Zhang } 1001c8db22f6SHong Zhang 10026973769bSHong Zhang #if defined(USE_ARRAY) 10036973769bSHong Zhang /* allocate an array for implementing sparse inner-product */ 1004e2cac8caSJed Brown ierr = PetscMalloc((A->cmap->n+1)*sizeof(MatScalar),&spdot);CHKERRQ(ierr); 1005e2cac8caSJed Brown ierr = PetscMemzero(spdot,(A->cmap->n+1)*sizeof(MatScalar));CHKERRQ(ierr); 10066973769bSHong Zhang #endif 10076973769bSHong Zhang 100881d82fe4SHong Zhang /* clear old values in C */ 1009aa1aec99SHong Zhang if (!c->a){ 1010aa1aec99SHong Zhang ierr = PetscMalloc((ci[cm]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr); 1011aa1aec99SHong Zhang c->a = ca; 1012aa1aec99SHong Zhang c->free_a = PETSC_TRUE; 1013aa1aec99SHong Zhang } else { 1014aa1aec99SHong Zhang ca = c->a; 1015aa1aec99SHong Zhang } 101681d82fe4SHong Zhang ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr); 10175df89d91SHong Zhang 10181fa1209cSHong Zhang for (i=0; i<cm; i++) { 101981d82fe4SHong Zhang anzi = ai[i+1] - ai[i]; 102081d82fe4SHong Zhang acol = aj + ai[i]; 10216973769bSHong Zhang aval = aa + ai[i]; 10221fa1209cSHong Zhang cnzi = ci[i+1] - ci[i]; 10231fa1209cSHong Zhang ccol = cj + ci[i]; 10246973769bSHong Zhang cval = ca + ci[i]; 10251fa1209cSHong Zhang for (j=0; j<cnzi; j++){ 102681d82fe4SHong Zhang brow = ccol[j]; 102781d82fe4SHong Zhang bnzj = bi[brow+1] - bi[brow]; 102881d82fe4SHong Zhang bcol = bj + bi[brow]; 10296973769bSHong Zhang bval = ba + bi[brow]; 10306973769bSHong Zhang 103181d82fe4SHong Zhang /* perform sparse inner-product c(i,j)=A[i,:]*B[j,:]^T */ 10326973769bSHong Zhang #if defined(USE_ARRAY) 10336973769bSHong Zhang /* put ba to spdot array */ 10346973769bSHong Zhang for (nextb=0; nextb<bnzj; nextb++) spdot[bcol[nextb]] = bval[nextb]; 10356973769bSHong Zhang /* c(i,j)=A[i,:]*B[j,:]^T */ 10366973769bSHong Zhang for (nexta=0; nexta<anzi; nexta++){ 10376973769bSHong Zhang cval[j] += spdot[acol[nexta]]*aval[nexta]; 10386973769bSHong Zhang } 10396973769bSHong Zhang /* zero spdot array */ 10406973769bSHong Zhang for (nextb=0; nextb<bnzj; nextb++) spdot[bcol[nextb]] = 0.0; 10416973769bSHong Zhang #else 104281d82fe4SHong Zhang nexta = 0; nextb = 0; 104381d82fe4SHong Zhang while (nexta<anzi && nextb<bnzj){ 104481d82fe4SHong Zhang while (acol[nexta] < bcol[nextb] && nexta < anzi) nexta++; 104581d82fe4SHong Zhang if (nexta == anzi) break; 104681d82fe4SHong Zhang while (acol[nexta] > bcol[nextb] && nextb < bnzj) nextb++; 104781d82fe4SHong Zhang if (nextb == bnzj) break; 104881d82fe4SHong Zhang if (acol[nexta] == bcol[nextb]){ 10496973769bSHong Zhang cval[j] += aval[nexta]*bval[nextb]; 105081d82fe4SHong Zhang nexta++; nextb++; 105181d82fe4SHong Zhang flops += 2; 10521fa1209cSHong Zhang } 10531fa1209cSHong Zhang } 10546973769bSHong Zhang #endif 105581d82fe4SHong Zhang } 105681d82fe4SHong Zhang } 105781d82fe4SHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 105881d82fe4SHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 105981d82fe4SHong Zhang ierr = PetscLogFlops(flops);CHKERRQ(ierr); 10606973769bSHong Zhang #if defined(USE_ARRAY) 10616973769bSHong Zhang ierr = PetscFree(spdot); 10626973769bSHong Zhang #endif 10635df89d91SHong Zhang PetscFunctionReturn(0); 10645df89d91SHong Zhang } 10655df89d91SHong Zhang 10665df89d91SHong Zhang #undef __FUNCT__ 106775648e8dSHong Zhang #define __FUNCT__ "MatTransposeMatMult_SeqAIJ_SeqAIJ" 106875648e8dSHong Zhang PetscErrorCode MatTransposeMatMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) { 10695df89d91SHong Zhang PetscErrorCode ierr; 10705df89d91SHong Zhang 10715df89d91SHong Zhang PetscFunctionBegin; 10725df89d91SHong Zhang if (scall == MAT_INITIAL_MATRIX){ 107375648e8dSHong Zhang ierr = MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr); 10745df89d91SHong Zhang } 107575648e8dSHong Zhang ierr = MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr); 10765df89d91SHong Zhang PetscFunctionReturn(0); 10775df89d91SHong Zhang } 10785df89d91SHong Zhang 10795df89d91SHong Zhang #undef __FUNCT__ 108075648e8dSHong Zhang #define __FUNCT__ "MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ" 108175648e8dSHong Zhang PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C) 10825df89d91SHong Zhang { 1083bc011b1eSHong Zhang PetscErrorCode ierr; 1084bc011b1eSHong Zhang Mat At; 108538baddfdSBarry Smith PetscInt *ati,*atj; 1086bc011b1eSHong Zhang 1087bc011b1eSHong Zhang PetscFunctionBegin; 1088bc011b1eSHong Zhang /* create symbolic At */ 1089bc011b1eSHong Zhang ierr = MatGetSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr); 1090d0f46423SBarry Smith ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,A->cmap->n,A->rmap->n,ati,atj,PETSC_NULL,&At);CHKERRQ(ierr); 1091a2f3521dSMark F. Adams At->rmap->bs = A->cmap->bs; 1092a2f3521dSMark F. Adams At->cmap->bs = B->cmap->bs; 1093bc011b1eSHong Zhang 1094bc011b1eSHong Zhang /* get symbolic C=At*B */ 1095bc011b1eSHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(At,B,fill,C);CHKERRQ(ierr); 1096bc011b1eSHong Zhang 1097bc011b1eSHong Zhang /* clean up */ 10986bf464f9SBarry Smith ierr = MatDestroy(&At);CHKERRQ(ierr); 1099bc011b1eSHong Zhang ierr = MatRestoreSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr); 1100bc011b1eSHong Zhang PetscFunctionReturn(0); 1101bc011b1eSHong Zhang } 1102bc011b1eSHong Zhang 1103bc011b1eSHong Zhang #undef __FUNCT__ 110475648e8dSHong Zhang #define __FUNCT__ "MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ" 110575648e8dSHong Zhang PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C) 1106bc011b1eSHong Zhang { 1107bc011b1eSHong Zhang PetscErrorCode ierr; 11080fbc74f4SHong Zhang Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data; 1109d0f46423SBarry Smith PetscInt am=A->rmap->n,anzi,*ai=a->i,*aj=a->j,*bi=b->i,*bj,bnzi,nextb; 1110d0f46423SBarry Smith PetscInt cm=C->rmap->n,*ci=c->i,*cj=c->j,crow,*cjj,i,j,k; 1111d13dce4bSSatish Balay PetscLogDouble flops=0.0; 1112aa1aec99SHong Zhang MatScalar *aa=a->a,*ba,*ca,*caj; 1113bc011b1eSHong Zhang 1114bc011b1eSHong Zhang PetscFunctionBegin; 1115aa1aec99SHong Zhang if (!c->a){ 1116aa1aec99SHong Zhang ierr = PetscMalloc((ci[cm]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr); 1117aa1aec99SHong Zhang c->a = ca; 1118aa1aec99SHong Zhang c->free_a = PETSC_TRUE; 1119aa1aec99SHong Zhang } else { 1120aa1aec99SHong Zhang ca = c->a; 1121aa1aec99SHong Zhang } 1122bc011b1eSHong Zhang /* clear old values in C */ 1123bc011b1eSHong Zhang ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr); 1124bc011b1eSHong Zhang 1125bc011b1eSHong Zhang /* compute A^T*B using outer product (A^T)[:,i]*B[i,:] */ 1126bc011b1eSHong Zhang for (i=0;i<am;i++) { 1127bc011b1eSHong Zhang bj = b->j + bi[i]; 1128bc011b1eSHong Zhang ba = b->a + bi[i]; 1129bc011b1eSHong Zhang bnzi = bi[i+1] - bi[i]; 1130bc011b1eSHong Zhang anzi = ai[i+1] - ai[i]; 1131bc011b1eSHong Zhang for (j=0; j<anzi; j++) { 1132bc011b1eSHong Zhang nextb = 0; 11330fbc74f4SHong Zhang crow = *aj++; 1134bc011b1eSHong Zhang cjj = cj + ci[crow]; 1135bc011b1eSHong Zhang caj = ca + ci[crow]; 1136bc011b1eSHong Zhang /* perform sparse axpy operation. Note cjj includes bj. */ 1137bc011b1eSHong Zhang for (k=0; nextb<bnzi; k++) { 11380fbc74f4SHong Zhang if (cjj[k] == *(bj+nextb)) { /* ccol == bcol */ 11390fbc74f4SHong Zhang caj[k] += (*aa)*(*(ba+nextb)); 1140bc011b1eSHong Zhang nextb++; 1141bc011b1eSHong Zhang } 1142bc011b1eSHong Zhang } 1143bc011b1eSHong Zhang flops += 2*bnzi; 11440fbc74f4SHong Zhang aa++; 1145bc011b1eSHong Zhang } 1146bc011b1eSHong Zhang } 1147bc011b1eSHong Zhang 1148bc011b1eSHong Zhang /* Assemble the final matrix and clean up */ 1149bc011b1eSHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1150bc011b1eSHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1151bc011b1eSHong Zhang ierr = PetscLogFlops(flops);CHKERRQ(ierr); 1152bc011b1eSHong Zhang PetscFunctionReturn(0); 1153bc011b1eSHong Zhang } 11549123193aSHong Zhang 1155ec01deb9SMatthew Knepley EXTERN_C_BEGIN 11569123193aSHong Zhang #undef __FUNCT__ 11579123193aSHong Zhang #define __FUNCT__ "MatMatMult_SeqAIJ_SeqDense" 11589123193aSHong Zhang PetscErrorCode MatMatMult_SeqAIJ_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 11599123193aSHong Zhang { 11609123193aSHong Zhang PetscErrorCode ierr; 11619123193aSHong Zhang 11629123193aSHong Zhang PetscFunctionBegin; 11639123193aSHong Zhang if (scall == MAT_INITIAL_MATRIX){ 11649123193aSHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqDense(A,B,fill,C);CHKERRQ(ierr); 11659123193aSHong Zhang } 11669123193aSHong Zhang ierr = MatMatMultNumeric_SeqAIJ_SeqDense(A,B,*C);CHKERRQ(ierr); 11679123193aSHong Zhang PetscFunctionReturn(0); 11689123193aSHong Zhang } 1169ec01deb9SMatthew Knepley EXTERN_C_END 1170edd81eecSMatthew Knepley 11719123193aSHong Zhang #undef __FUNCT__ 11729123193aSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqDense" 11739123193aSHong Zhang PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C) 11749123193aSHong Zhang { 11759123193aSHong Zhang PetscErrorCode ierr; 11769123193aSHong Zhang 11779123193aSHong Zhang PetscFunctionBegin; 11785a586d82SBarry Smith ierr = MatMatMultSymbolic_SeqDense_SeqDense(A,B,0.0,C);CHKERRQ(ierr); 11798cdbd757SHong Zhang (*C)->ops->matmult = MatMatMult_SeqAIJ_SeqDense; 11809123193aSHong Zhang PetscFunctionReturn(0); 11819123193aSHong Zhang } 11829123193aSHong Zhang 11839123193aSHong Zhang #undef __FUNCT__ 11849123193aSHong Zhang #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqDense" 11859123193aSHong Zhang PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqDense(Mat A,Mat B,Mat C) 11869123193aSHong Zhang { 1187f32d5d43SBarry Smith Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data; 11889123193aSHong Zhang PetscErrorCode ierr; 1189dd6ea824SBarry Smith PetscScalar *b,*c,r1,r2,r3,r4,*b1,*b2,*b3,*b4; 1190dd6ea824SBarry Smith MatScalar *aa; 1191d0f46423SBarry Smith PetscInt cm=C->rmap->n, cn=B->cmap->n, bm=B->rmap->n, col, i,j,n,*aj, am = A->rmap->n; 1192f32d5d43SBarry Smith PetscInt am2 = 2*am, am3 = 3*am, bm4 = 4*bm,colam; 11939123193aSHong Zhang 11949123193aSHong Zhang PetscFunctionBegin; 1195f32d5d43SBarry Smith if (!cm || !cn) PetscFunctionReturn(0); 1196e32f2f54SBarry 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); 1197e32f2f54SBarry 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); 1198e32f2f54SBarry 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); 11998c778c55SBarry Smith ierr = MatDenseGetArray(B,&b);CHKERRQ(ierr); 12008c778c55SBarry Smith ierr = MatDenseGetArray(C,&c);CHKERRQ(ierr); 1201f32d5d43SBarry Smith b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm; 1202f32d5d43SBarry Smith for (col=0; col<cn-4; col += 4){ /* over columns of C */ 1203f32d5d43SBarry Smith colam = col*am; 1204f32d5d43SBarry Smith for (i=0; i<am; i++) { /* over rows of C in those columns */ 1205f32d5d43SBarry Smith r1 = r2 = r3 = r4 = 0.0; 1206f32d5d43SBarry Smith n = a->i[i+1] - a->i[i]; 1207f32d5d43SBarry Smith aj = a->j + a->i[i]; 1208f32d5d43SBarry Smith aa = a->a + a->i[i]; 1209f32d5d43SBarry Smith for (j=0; j<n; j++) { 1210f32d5d43SBarry Smith r1 += (*aa)*b1[*aj]; 1211f32d5d43SBarry Smith r2 += (*aa)*b2[*aj]; 1212f32d5d43SBarry Smith r3 += (*aa)*b3[*aj]; 1213f32d5d43SBarry Smith r4 += (*aa++)*b4[*aj++]; 12149123193aSHong Zhang } 1215f32d5d43SBarry Smith c[colam + i] = r1; 1216f32d5d43SBarry Smith c[colam + am + i] = r2; 1217f32d5d43SBarry Smith c[colam + am2 + i] = r3; 1218f32d5d43SBarry Smith c[colam + am3 + i] = r4; 1219f32d5d43SBarry Smith } 1220f32d5d43SBarry Smith b1 += bm4; 1221f32d5d43SBarry Smith b2 += bm4; 1222f32d5d43SBarry Smith b3 += bm4; 1223f32d5d43SBarry Smith b4 += bm4; 1224f32d5d43SBarry Smith } 1225f32d5d43SBarry Smith for (;col<cn; col++){ /* over extra columns of C */ 1226f32d5d43SBarry Smith for (i=0; i<am; i++) { /* over rows of C in those columns */ 1227f32d5d43SBarry Smith r1 = 0.0; 1228f32d5d43SBarry Smith n = a->i[i+1] - a->i[i]; 1229f32d5d43SBarry Smith aj = a->j + a->i[i]; 1230f32d5d43SBarry Smith aa = a->a + a->i[i]; 1231f32d5d43SBarry Smith 1232f32d5d43SBarry Smith for (j=0; j<n; j++) { 1233f32d5d43SBarry Smith r1 += (*aa++)*b1[*aj++]; 1234f32d5d43SBarry Smith } 1235f32d5d43SBarry Smith c[col*am + i] = r1; 1236f32d5d43SBarry Smith } 1237f32d5d43SBarry Smith b1 += bm; 1238f32d5d43SBarry Smith } 1239dc0b31edSSatish Balay ierr = PetscLogFlops(cn*(2.0*a->nz));CHKERRQ(ierr); 12408c778c55SBarry Smith ierr = MatDenseRestoreArray(B,&b);CHKERRQ(ierr); 12418c778c55SBarry Smith ierr = MatDenseRestoreArray(C,&c);CHKERRQ(ierr); 1242f32d5d43SBarry Smith ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1243f32d5d43SBarry Smith ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1244f32d5d43SBarry Smith PetscFunctionReturn(0); 1245f32d5d43SBarry Smith } 1246f32d5d43SBarry Smith 1247f32d5d43SBarry Smith /* 12484324174eSBarry Smith Note very similar to MatMult_SeqAIJ(), should generate both codes from same base 1249f32d5d43SBarry Smith */ 1250f32d5d43SBarry Smith #undef __FUNCT__ 1251f32d5d43SBarry Smith #define __FUNCT__ "MatMatMultNumericAdd_SeqAIJ_SeqDense" 1252f32d5d43SBarry Smith PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat A,Mat B,Mat C) 1253f32d5d43SBarry Smith { 1254f32d5d43SBarry Smith Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data; 1255f32d5d43SBarry Smith PetscErrorCode ierr; 1256dd6ea824SBarry Smith PetscScalar *b,*c,r1,r2,r3,r4,*b1,*b2,*b3,*b4; 1257dd6ea824SBarry Smith MatScalar *aa; 1258d0f46423SBarry 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; 12594324174eSBarry Smith PetscInt am2 = 2*am, am3 = 3*am, bm4 = 4*bm,colam,*ridx; 1260f32d5d43SBarry Smith 1261f32d5d43SBarry Smith PetscFunctionBegin; 1262f32d5d43SBarry Smith if (!cm || !cn) PetscFunctionReturn(0); 12638c778c55SBarry Smith ierr = MatDenseGetArray(B,&b);CHKERRQ(ierr); 12648c778c55SBarry Smith ierr = MatDenseGetArray(C,&c);CHKERRQ(ierr); 1265f32d5d43SBarry Smith b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm; 12664324174eSBarry Smith 12674324174eSBarry Smith if (a->compressedrow.use){ /* use compressed row format */ 12684324174eSBarry Smith for (col=0; col<cn-4; col += 4){ /* over columns of C */ 12694324174eSBarry Smith colam = col*am; 12704324174eSBarry Smith arm = a->compressedrow.nrows; 12714324174eSBarry Smith ii = a->compressedrow.i; 12724324174eSBarry Smith ridx = a->compressedrow.rindex; 12734324174eSBarry Smith for (i=0; i<arm; i++) { /* over rows of C in those columns */ 12744324174eSBarry Smith r1 = r2 = r3 = r4 = 0.0; 12754324174eSBarry Smith n = ii[i+1] - ii[i]; 12764324174eSBarry Smith aj = a->j + ii[i]; 12774324174eSBarry Smith aa = a->a + ii[i]; 12784324174eSBarry Smith for (j=0; j<n; j++) { 12794324174eSBarry Smith r1 += (*aa)*b1[*aj]; 12804324174eSBarry Smith r2 += (*aa)*b2[*aj]; 12814324174eSBarry Smith r3 += (*aa)*b3[*aj]; 12824324174eSBarry Smith r4 += (*aa++)*b4[*aj++]; 12834324174eSBarry Smith } 12844324174eSBarry Smith c[colam + ridx[i]] += r1; 12854324174eSBarry Smith c[colam + am + ridx[i]] += r2; 12864324174eSBarry Smith c[colam + am2 + ridx[i]] += r3; 12874324174eSBarry Smith c[colam + am3 + ridx[i]] += r4; 12884324174eSBarry Smith } 12894324174eSBarry Smith b1 += bm4; 12904324174eSBarry Smith b2 += bm4; 12914324174eSBarry Smith b3 += bm4; 12924324174eSBarry Smith b4 += bm4; 12934324174eSBarry Smith } 12944324174eSBarry Smith for (;col<cn; col++){ /* over extra columns of C */ 12954324174eSBarry Smith colam = col*am; 12964324174eSBarry Smith arm = a->compressedrow.nrows; 12974324174eSBarry Smith ii = a->compressedrow.i; 12984324174eSBarry Smith ridx = a->compressedrow.rindex; 12994324174eSBarry Smith for (i=0; i<arm; i++) { /* over rows of C in those columns */ 13004324174eSBarry Smith r1 = 0.0; 13014324174eSBarry Smith n = ii[i+1] - ii[i]; 13024324174eSBarry Smith aj = a->j + ii[i]; 13034324174eSBarry Smith aa = a->a + ii[i]; 13044324174eSBarry Smith 13054324174eSBarry Smith for (j=0; j<n; j++) { 13064324174eSBarry Smith r1 += (*aa++)*b1[*aj++]; 13074324174eSBarry Smith } 13084324174eSBarry Smith c[col*am + ridx[i]] += r1; 13094324174eSBarry Smith } 13104324174eSBarry Smith b1 += bm; 13114324174eSBarry Smith } 13124324174eSBarry Smith } else { 1313f32d5d43SBarry Smith for (col=0; col<cn-4; col += 4){ /* over columns of C */ 1314f32d5d43SBarry Smith colam = col*am; 1315f32d5d43SBarry Smith for (i=0; i<am; i++) { /* over rows of C in those columns */ 1316f32d5d43SBarry Smith r1 = r2 = r3 = r4 = 0.0; 1317f32d5d43SBarry Smith n = a->i[i+1] - a->i[i]; 1318f32d5d43SBarry Smith aj = a->j + a->i[i]; 1319f32d5d43SBarry Smith aa = a->a + a->i[i]; 1320f32d5d43SBarry Smith for (j=0; j<n; j++) { 1321f32d5d43SBarry Smith r1 += (*aa)*b1[*aj]; 1322f32d5d43SBarry Smith r2 += (*aa)*b2[*aj]; 1323f32d5d43SBarry Smith r3 += (*aa)*b3[*aj]; 1324f32d5d43SBarry Smith r4 += (*aa++)*b4[*aj++]; 1325f32d5d43SBarry Smith } 1326f32d5d43SBarry Smith c[colam + i] += r1; 1327f32d5d43SBarry Smith c[colam + am + i] += r2; 1328f32d5d43SBarry Smith c[colam + am2 + i] += r3; 1329f32d5d43SBarry Smith c[colam + am3 + i] += r4; 1330f32d5d43SBarry Smith } 1331f32d5d43SBarry Smith b1 += bm4; 1332f32d5d43SBarry Smith b2 += bm4; 1333f32d5d43SBarry Smith b3 += bm4; 1334f32d5d43SBarry Smith b4 += bm4; 1335f32d5d43SBarry Smith } 1336f32d5d43SBarry Smith for (;col<cn; col++){ /* over extra columns of C */ 1337f32d5d43SBarry Smith for (i=0; i<am; i++) { /* over rows of C in those columns */ 1338f32d5d43SBarry Smith r1 = 0.0; 1339f32d5d43SBarry Smith n = a->i[i+1] - a->i[i]; 1340f32d5d43SBarry Smith aj = a->j + a->i[i]; 1341f32d5d43SBarry Smith aa = a->a + a->i[i]; 1342f32d5d43SBarry Smith 1343f32d5d43SBarry Smith for (j=0; j<n; j++) { 1344f32d5d43SBarry Smith r1 += (*aa++)*b1[*aj++]; 1345f32d5d43SBarry Smith } 1346f32d5d43SBarry Smith c[col*am + i] += r1; 1347f32d5d43SBarry Smith } 1348f32d5d43SBarry Smith b1 += bm; 1349f32d5d43SBarry Smith } 13504324174eSBarry Smith } 1351dc0b31edSSatish Balay ierr = PetscLogFlops(cn*2.0*a->nz);CHKERRQ(ierr); 13528c778c55SBarry Smith ierr = MatDenseRestoreArray(B,&b);CHKERRQ(ierr); 13538c778c55SBarry Smith ierr = MatDenseRestoreArray(C,&c);CHKERRQ(ierr); 13549123193aSHong Zhang PetscFunctionReturn(0); 13559123193aSHong Zhang } 1356b1683b59SHong Zhang 1357b1683b59SHong Zhang #undef __FUNCT__ 1358b9af6bddSHong Zhang #define __FUNCT__ "MatTransColoringApplySpToDen_SeqAIJ" 1359b9af6bddSHong Zhang PetscErrorCode MatTransColoringApplySpToDen_SeqAIJ(MatTransposeColoring coloring,Mat B,Mat Btdense) 1360c8db22f6SHong Zhang { 1361c8db22f6SHong Zhang PetscErrorCode ierr; 13622f5f1f90SHong Zhang Mat_SeqAIJ *b = (Mat_SeqAIJ*)B->data; 13632f5f1f90SHong Zhang Mat_SeqDense *btdense = (Mat_SeqDense*)Btdense->data; 13642f5f1f90SHong Zhang PetscInt *bi=b->i,*bj=b->j; 13652f5f1f90SHong Zhang PetscInt m=Btdense->rmap->n,n=Btdense->cmap->n,j,k,l,col,anz,*btcol,brow,ncolumns; 13662f5f1f90SHong Zhang MatScalar *btval,*btval_den,*ba=b->a; 13672f5f1f90SHong Zhang PetscInt *columns=coloring->columns,*colorforcol=coloring->colorforcol,ncolors=coloring->ncolors; 1368c8db22f6SHong Zhang 1369c8db22f6SHong Zhang PetscFunctionBegin; 13702f5f1f90SHong Zhang btval_den=btdense->v; 13712f5f1f90SHong Zhang ierr = PetscMemzero(btval_den,(m*n)*sizeof(MatScalar));CHKERRQ(ierr); 13722f5f1f90SHong Zhang for (k=0; k<ncolors; k++) { 13732f5f1f90SHong Zhang ncolumns = coloring->ncolumns[k]; 13742f5f1f90SHong Zhang for (l=0; l<ncolumns; l++) { /* insert a row of B to a column of Btdense */ 1375d2853540SHong Zhang col = *(columns + colorforcol[k] + l); 13762f5f1f90SHong Zhang btcol = bj + bi[col]; 13772f5f1f90SHong Zhang btval = ba + bi[col]; 13782f5f1f90SHong Zhang anz = bi[col+1] - bi[col]; 1379d2853540SHong Zhang for (j=0; j<anz; j++){ 13802f5f1f90SHong Zhang brow = btcol[j]; 13812f5f1f90SHong Zhang btval_den[brow] = btval[j]; 1382c8db22f6SHong Zhang } 1383c8db22f6SHong Zhang } 13842f5f1f90SHong Zhang btval_den += m; 1385c8db22f6SHong Zhang } 1386c8db22f6SHong Zhang PetscFunctionReturn(0); 1387c8db22f6SHong Zhang } 1388c8db22f6SHong Zhang 1389c8db22f6SHong Zhang #undef __FUNCT__ 1390b9af6bddSHong Zhang #define __FUNCT__ "MatTransColoringApplyDenToSp_SeqAIJ" 1391b9af6bddSHong Zhang PetscErrorCode MatTransColoringApplyDenToSp_SeqAIJ(MatTransposeColoring matcoloring,Mat Cden,Mat Csp) 13928972f759SHong Zhang { 13938972f759SHong Zhang PetscErrorCode ierr; 1394b2d2b04fSHong Zhang Mat_SeqAIJ *csp = (Mat_SeqAIJ*)Csp->data; 13952f5f1f90SHong Zhang PetscInt k,l,*row,*idx,m,ncolors=matcoloring->ncolors,nrows; 1396b2d2b04fSHong Zhang PetscScalar *ca_den,*cp_den,*ca=csp->a; 13972f5f1f90SHong Zhang PetscInt *rows=matcoloring->rows,*spidx=matcoloring->columnsforspidx,*colorforrow=matcoloring->colorforrow; 13988972f759SHong Zhang 13998972f759SHong Zhang PetscFunctionBegin; 14008972f759SHong Zhang ierr = MatGetLocalSize(Csp,&m,PETSC_NULL);CHKERRQ(ierr); 1401*a3fe58edSHong Zhang ierr = MatDenseGetArray(Cden,&ca_den);CHKERRQ(ierr); 1402b2d2b04fSHong Zhang cp_den = ca_den; 14032f5f1f90SHong Zhang for (k=0; k<ncolors; k++) { 14042f5f1f90SHong Zhang nrows = matcoloring->nrows[k]; 14052f5f1f90SHong Zhang row = rows + colorforrow[k]; 14062f5f1f90SHong Zhang idx = spidx + colorforrow[k]; 14072f5f1f90SHong Zhang for (l=0; l<nrows; l++){ 14082f5f1f90SHong Zhang ca[idx[l]] = cp_den[row[l]]; 1409b2d2b04fSHong Zhang } 1410b2d2b04fSHong Zhang cp_den += m; 1411b2d2b04fSHong Zhang } 1412*a3fe58edSHong Zhang ierr = MatDenseRestoreArray(Cden,&ca_den);CHKERRQ(ierr); 14138972f759SHong Zhang PetscFunctionReturn(0); 14148972f759SHong Zhang } 14158972f759SHong Zhang 1416e99be685SHong Zhang /* 1417e99be685SHong Zhang MatGetColumnIJ_SeqAIJ_Color() and MatRestoreColumnIJ_SeqAIJ_Color() are customized from 1418e99be685SHong Zhang MatGetColumnIJ_SeqAIJ() and MatRestoreColumnIJ_SeqAIJ() by adding an output 1419e99be685SHong Zhang spidx[], index of a->j, to be used for setting 'columnsforspidx' in MatTransposeColoringCreate_SeqAIJ(). 1420e99be685SHong Zhang */ 1421e99be685SHong Zhang #undef __FUNCT__ 1422e99be685SHong Zhang #define __FUNCT__ "MatGetColumnIJ_SeqAIJ_Color" 1423e99be685SHong Zhang PetscErrorCode MatGetColumnIJ_SeqAIJ_Color(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool inodecompressed,PetscInt *nn,PetscInt *ia[],PetscInt *ja[],PetscInt *spidx[],PetscBool *done) 1424e99be685SHong Zhang { 1425e99be685SHong Zhang Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 1426e99be685SHong Zhang PetscErrorCode ierr; 1427e99be685SHong Zhang PetscInt i,*collengths,*cia,*cja,n = A->cmap->n,m = A->rmap->n; 1428e99be685SHong Zhang PetscInt nz = a->i[m],row,*jj,mr,col; 1429e99be685SHong Zhang PetscInt *cspidx; 1430e99be685SHong Zhang 1431e99be685SHong Zhang PetscFunctionBegin; 1432e99be685SHong Zhang *nn = n; 1433e99be685SHong Zhang if (!ia) PetscFunctionReturn(0); 1434e99be685SHong Zhang if (symmetric) { 1435e99be685SHong Zhang SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"MatGetColumnIJ_SeqAIJ_Color() not supported for the case symmetric"); 1436e99be685SHong Zhang ierr = MatToSymmetricIJ_SeqAIJ(A->rmap->n,a->i,a->j,0,oshift,ia,ja);CHKERRQ(ierr); 1437e99be685SHong Zhang } else { 1438e99be685SHong Zhang ierr = PetscMalloc((n+1)*sizeof(PetscInt),&collengths);CHKERRQ(ierr); 1439e99be685SHong Zhang ierr = PetscMemzero(collengths,n*sizeof(PetscInt));CHKERRQ(ierr); 1440e99be685SHong Zhang ierr = PetscMalloc((n+1)*sizeof(PetscInt),&cia);CHKERRQ(ierr); 1441e99be685SHong Zhang ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&cja);CHKERRQ(ierr); 1442e99be685SHong Zhang ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&cspidx);CHKERRQ(ierr); 1443e99be685SHong Zhang jj = a->j; 1444e99be685SHong Zhang for (i=0; i<nz; i++) { 1445e99be685SHong Zhang collengths[jj[i]]++; 1446e99be685SHong Zhang } 1447e99be685SHong Zhang cia[0] = oshift; 1448e99be685SHong Zhang for (i=0; i<n; i++) { 1449e99be685SHong Zhang cia[i+1] = cia[i] + collengths[i]; 1450e99be685SHong Zhang } 1451e99be685SHong Zhang ierr = PetscMemzero(collengths,n*sizeof(PetscInt));CHKERRQ(ierr); 1452e99be685SHong Zhang jj = a->j; 1453e99be685SHong Zhang for (row=0; row<m; row++) { 1454e99be685SHong Zhang mr = a->i[row+1] - a->i[row]; 1455e99be685SHong Zhang for (i=0; i<mr; i++) { 1456e99be685SHong Zhang col = *jj++; 1457e99be685SHong Zhang cspidx[cia[col] + collengths[col] - oshift] = a->i[row] + i; /* index of a->j */ 1458e99be685SHong Zhang cja[cia[col] + collengths[col]++ - oshift] = row + oshift; 1459e99be685SHong Zhang } 1460e99be685SHong Zhang } 1461e99be685SHong Zhang ierr = PetscFree(collengths);CHKERRQ(ierr); 1462e99be685SHong Zhang *ia = cia; *ja = cja; 1463e99be685SHong Zhang *spidx = cspidx; 1464e99be685SHong Zhang } 1465e99be685SHong Zhang PetscFunctionReturn(0); 1466e99be685SHong Zhang } 1467e99be685SHong Zhang 1468e99be685SHong Zhang #undef __FUNCT__ 1469e99be685SHong Zhang #define __FUNCT__ "MatRestoreColumnIJ_SeqAIJ_Color" 1470e99be685SHong Zhang PetscErrorCode MatRestoreColumnIJ_SeqAIJ_Color(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool inodecompressed,PetscInt *n,PetscInt *ia[],PetscInt *ja[],PetscInt *spidx[],PetscBool *done) 1471e99be685SHong Zhang { 1472e99be685SHong Zhang PetscErrorCode ierr; 1473e99be685SHong Zhang 1474e99be685SHong Zhang PetscFunctionBegin; 1475e99be685SHong Zhang if (!ia) PetscFunctionReturn(0); 1476e99be685SHong Zhang 1477e99be685SHong Zhang ierr = PetscFree(*ia);CHKERRQ(ierr); 1478e99be685SHong Zhang ierr = PetscFree(*ja);CHKERRQ(ierr); 1479e99be685SHong Zhang ierr = PetscFree(*spidx);CHKERRQ(ierr); 1480e99be685SHong Zhang PetscFunctionReturn(0); 1481e99be685SHong Zhang } 1482e99be685SHong Zhang 14838972f759SHong Zhang #undef __FUNCT__ 1484b9af6bddSHong Zhang #define __FUNCT__ "MatTransposeColoringCreate_SeqAIJ" 1485b9af6bddSHong Zhang PetscErrorCode MatTransposeColoringCreate_SeqAIJ(Mat mat,ISColoring iscoloring,MatTransposeColoring c) 1486b1683b59SHong Zhang { 1487b1683b59SHong Zhang PetscErrorCode ierr; 1488d2853540SHong Zhang PetscInt i,n,nrows,N,j,k,m,*row_idx,*ci,*cj,ncols,col,cm; 1489b1683b59SHong Zhang const PetscInt *is; 1490b28d80bdSHong Zhang PetscInt nis = iscoloring->n,*rowhit,bs = 1; 1491b1683b59SHong Zhang IS *isa; 1492b1683b59SHong Zhang PetscBool flg1,flg2; 1493b28d80bdSHong Zhang Mat_SeqAIJ *csp = (Mat_SeqAIJ*)mat->data; 1494e99be685SHong Zhang PetscInt *colorforrow,*rows,*rows_i,*columnsforspidx,*columnsforspidx_i,*idxhit,*spidx; 1495d2853540SHong Zhang PetscInt *colorforcol,*columns,*columns_i; 1496b1683b59SHong Zhang 1497b1683b59SHong Zhang PetscFunctionBegin; 1498b1683b59SHong Zhang ierr = ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);CHKERRQ(ierr); 1499e99be685SHong Zhang 1500b1683b59SHong Zhang /* this is ugly way to get blocksize but cannot call MatGetBlockSize() because AIJ can have bs > 1 */ 1501251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)mat,MATSEQBAIJ,&flg1);CHKERRQ(ierr); 1502251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)mat,MATMPIBAIJ,&flg2);CHKERRQ(ierr); 1503b1683b59SHong Zhang if (flg1 || flg2) { 1504b1683b59SHong Zhang ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr); 1505b1683b59SHong Zhang } 1506b1683b59SHong Zhang 1507b1683b59SHong Zhang N = mat->cmap->N/bs; 1508b1683b59SHong Zhang c->M = mat->rmap->N/bs; /* set total rows, columns and local rows */ 1509b1683b59SHong Zhang c->N = mat->cmap->N/bs; 1510b1683b59SHong Zhang c->m = mat->rmap->N/bs; 1511b1683b59SHong Zhang c->rstart = 0; 1512b1683b59SHong Zhang 1513b1683b59SHong Zhang c->ncolors = nis; 1514b1683b59SHong Zhang ierr = PetscMalloc(nis*sizeof(PetscInt),&c->ncolumns);CHKERRQ(ierr); 1515b28d80bdSHong Zhang ierr = PetscMalloc(nis*sizeof(PetscInt),&c->nrows);CHKERRQ(ierr); 1516d2853540SHong Zhang ierr = PetscMalloc2(csp->nz+1,PetscInt,&rows,csp->nz+1,PetscInt,&columnsforspidx);CHKERRQ(ierr); 1517d2853540SHong Zhang ierr = PetscMalloc((nis+1)*sizeof(PetscInt),&colorforrow);CHKERRQ(ierr); 1518d2853540SHong Zhang colorforrow[0] = 0; 1519d2853540SHong Zhang rows_i = rows; 1520d2853540SHong Zhang columnsforspidx_i = columnsforspidx; 1521b1683b59SHong Zhang 1522d2853540SHong Zhang ierr = PetscMalloc((nis+1)*sizeof(PetscInt),&colorforcol);CHKERRQ(ierr); 1523d2853540SHong Zhang ierr = PetscMalloc((N+1)*sizeof(PetscInt),&columns);CHKERRQ(ierr); 1524d2853540SHong Zhang colorforcol[0] = 0; 1525d2853540SHong Zhang columns_i = columns; 1526d2853540SHong Zhang 1527f2691f8dSHong Zhang ierr = MatGetColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,PETSC_NULL);CHKERRQ(ierr); /* column-wise storage of mat */ 1528b1683b59SHong Zhang 1529b28d80bdSHong Zhang cm = c->m; 1530b28d80bdSHong Zhang ierr = PetscMalloc((cm+1)*sizeof(PetscInt),&rowhit);CHKERRQ(ierr); 1531e99be685SHong Zhang ierr = PetscMalloc((cm+1)*sizeof(PetscInt),&idxhit);CHKERRQ(ierr); 1532b1683b59SHong Zhang for (i=0; i<nis; i++) { 1533b1683b59SHong Zhang ierr = ISGetLocalSize(isa[i],&n);CHKERRQ(ierr); 1534b1683b59SHong Zhang ierr = ISGetIndices(isa[i],&is);CHKERRQ(ierr); 1535b1683b59SHong Zhang c->ncolumns[i] = n; 1536b1683b59SHong Zhang if (n) { 1537d2853540SHong Zhang ierr = PetscMemcpy(columns_i,is,n*sizeof(PetscInt));CHKERRQ(ierr); 1538b1683b59SHong Zhang } 1539d2853540SHong Zhang colorforcol[i+1] = colorforcol[i] + n; 1540d2853540SHong Zhang columns_i += n; 1541b1683b59SHong Zhang 1542b1683b59SHong Zhang /* fast, crude version requires O(N*N) work */ 1543e8354b3eSHong Zhang ierr = PetscMemzero(rowhit,cm*sizeof(PetscInt));CHKERRQ(ierr); 1544e99be685SHong Zhang 1545b1683b59SHong Zhang /* loop over columns*/ 1546b1683b59SHong Zhang for (j=0; j<n; j++) { 1547b1683b59SHong Zhang col = is[j]; 1548d2853540SHong Zhang row_idx = cj + ci[col]; 1549b1683b59SHong Zhang m = ci[col+1] - ci[col]; 1550b1683b59SHong Zhang /* loop over columns marking them in rowhit */ 1551b1683b59SHong Zhang for (k=0; k<m; k++) { 1552e99be685SHong Zhang idxhit[*row_idx] = spidx[ci[col] + k]; 1553d2853540SHong Zhang rowhit[*row_idx++] = col + 1; 1554b1683b59SHong Zhang } 1555b1683b59SHong Zhang } 1556b1683b59SHong Zhang /* count the number of hits */ 1557b1683b59SHong Zhang nrows = 0; 1558e8354b3eSHong Zhang for (j=0; j<cm; j++) { 1559b1683b59SHong Zhang if (rowhit[j]) nrows++; 1560b1683b59SHong Zhang } 1561b1683b59SHong Zhang c->nrows[i] = nrows; 1562d2853540SHong Zhang colorforrow[i+1] = colorforrow[i] + nrows; 1563d2853540SHong Zhang 1564b1683b59SHong Zhang nrows = 0; 1565e8354b3eSHong Zhang for (j=0; j<cm; j++) { 1566b1683b59SHong Zhang if (rowhit[j]) { 1567d2853540SHong Zhang rows_i[nrows] = j; 1568e99be685SHong Zhang columnsforspidx_i[nrows] = idxhit[j]; 1569b1683b59SHong Zhang nrows++; 1570b1683b59SHong Zhang } 1571b1683b59SHong Zhang } 1572b1683b59SHong Zhang ierr = ISRestoreIndices(isa[i],&is);CHKERRQ(ierr); 1573d2853540SHong Zhang rows_i += nrows; columnsforspidx_i += nrows; 1574b1683b59SHong Zhang } 1575f2691f8dSHong Zhang ierr = MatRestoreColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,PETSC_NULL);CHKERRQ(ierr); 1576b28d80bdSHong Zhang ierr = PetscFree(rowhit);CHKERRQ(ierr); 1577b1683b59SHong Zhang ierr = ISColoringRestoreIS(iscoloring,&isa);CHKERRQ(ierr); 1578d2853540SHong Zhang #if defined(PETSC_USE_DEBUG) 1579d2853540SHong Zhang if (csp->nz != colorforrow[nis]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"csp->nz %d != colorforrow[nis] %d",csp->nz,colorforrow[nis]); 1580d2853540SHong Zhang #endif 1581b28d80bdSHong Zhang 1582d2853540SHong Zhang c->colorforrow = colorforrow; 1583d2853540SHong Zhang c->rows = rows; 1584d2853540SHong Zhang c->columnsforspidx = columnsforspidx; 1585d2853540SHong Zhang c->colorforcol = colorforcol; 1586d2853540SHong Zhang c->columns = columns; 1587f94ccd6cSHong Zhang ierr = PetscFree(idxhit);CHKERRQ(ierr); 1588b1683b59SHong Zhang PetscFunctionReturn(0); 1589b1683b59SHong Zhang } 1590