1be1d678aSKris Buschelman 2d50806bdSBarry Smith /* 32499ec78SHong Zhang Defines matrix-matrix product routines for pairs of SeqAIJ matrices 4d50806bdSBarry Smith C = A * B 5d50806bdSBarry Smith */ 6d50806bdSBarry Smith 7c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/aij.h> /*I "petscmat.h" I*/ 8c6db04a5SJed Brown #include <../src/mat/utils/freespace.h> 9c6db04a5SJed Brown #include <petscbt.h> 10af0996ceSBarry Smith #include <petsc/private/isimpl.h> 1107475bc1SBarry Smith #include <../src/mat/impls/dense/seq/dense.h> 127bab7c10SHong Zhang 1358cf0668SJed Brown static PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_LLCondensed(Mat,Mat,PetscReal,Mat*); 14cd093f1dSJed Brown 155e5acdf2Sstefano_zampini #if defined(PETSC_HAVE_HYPRE) 165e5acdf2Sstefano_zampini PETSC_INTERN PetscErrorCode MatMatMultSymbolic_AIJ_AIJ_wHYPRE(Mat,Mat,PetscReal,Mat*); 175e5acdf2Sstefano_zampini #endif 185e5acdf2Sstefano_zampini 19150d2497SBarry Smith PETSC_INTERN PetscErrorCode MatMatMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 2038baddfdSBarry Smith { 21dfbe8321SBarry Smith PetscErrorCode ierr; 225e5acdf2Sstefano_zampini #if !defined(PETSC_HAVE_HYPRE) 23d013fc79Sandi selinger const char *algTypes[7] = {"sorted","scalable","scalable_fast","heap","btheap","llcondensed","combined"}; 245e5acdf2Sstefano_zampini PetscInt nalg = 7; 25d013fc79Sandi selinger #else 26d013fc79Sandi selinger const char *algTypes[8] = {"sorted","scalable","scalable_fast","heap","btheap","llcondensed","combined","hypre"}; 27d013fc79Sandi selinger PetscInt nalg = 8; 285e5acdf2Sstefano_zampini #endif 296540a9cdSHong Zhang PetscInt alg = 0; /* set default algorithm */ 30d013fc79Sandi selinger PetscBool combined = PETSC_FALSE; /* Indicates whether the symbolic stage already computed the numerical values. */ 315c66b693SKris Buschelman 325c66b693SKris Buschelman PetscFunctionBegin; 3326be0446SHong Zhang if (scall == MAT_INITIAL_MATRIX) { 34152983bfSHong Zhang ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr); 3568eacb73SHong Zhang PetscOptionsObject->alreadyprinted = PETSC_FALSE; /* a hack to ensure the option shows in '-help' */ 365e5acdf2Sstefano_zampini ierr = PetscOptionsEList("-matmatmult_via","Algorithmic approach","MatMatMult",algTypes,nalg,algTypes[0],&alg,NULL);CHKERRQ(ierr); 37d8bbc50fSBarry Smith ierr = PetscOptionsEnd();CHKERRQ(ierr); 383ff4c91cSHong Zhang ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr); 396540a9cdSHong Zhang switch (alg) { 406540a9cdSHong Zhang case 1: 41aacf854cSBarry Smith ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(A,B,fill,C);CHKERRQ(ierr); 426540a9cdSHong Zhang break; 436540a9cdSHong Zhang case 2: 446540a9cdSHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast(A,B,fill,C);CHKERRQ(ierr); 456540a9cdSHong Zhang break; 466540a9cdSHong Zhang case 3: 470ced3a2bSJed Brown ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_Heap(A,B,fill,C);CHKERRQ(ierr); 486540a9cdSHong Zhang break; 496540a9cdSHong Zhang case 4: 508a07c6f1SJed Brown ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_BTHeap(A,B,fill,C);CHKERRQ(ierr); 516540a9cdSHong Zhang break; 526540a9cdSHong Zhang case 5: 5358cf0668SJed Brown ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_LLCondensed(A,B,fill,C);CHKERRQ(ierr); 546540a9cdSHong Zhang break; 555e5acdf2Sstefano_zampini case 6: 56d013fc79Sandi selinger ierr = MatMatMult_SeqAIJ_SeqAIJ_Combined(A,B,fill,C);CHKERRQ(ierr); 57d013fc79Sandi selinger combined = PETSC_TRUE; 58d013fc79Sandi selinger break; 59d013fc79Sandi selinger #if defined(PETSC_HAVE_HYPRE) 60d013fc79Sandi selinger case 7: 615e5acdf2Sstefano_zampini ierr = MatMatMultSymbolic_AIJ_AIJ_wHYPRE(A,B,fill,C);CHKERRQ(ierr); 625e5acdf2Sstefano_zampini break; 635e5acdf2Sstefano_zampini #endif 646540a9cdSHong Zhang default: 6526be0446SHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr); 666540a9cdSHong Zhang break; 6725023028SHong Zhang } 683ff4c91cSHong Zhang ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr); 6926be0446SHong Zhang } 705c913ed7SHong Zhang 713ff4c91cSHong Zhang ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr); 72d013fc79Sandi selinger if (!combined) { 7301e47db4SHong Zhang ierr = (*(*C)->ops->matmultnumeric)(A,B,*C);CHKERRQ(ierr); 74d013fc79Sandi selinger } 753ff4c91cSHong Zhang ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr); 765c66b693SKris Buschelman PetscFunctionReturn(0); 775c66b693SKris Buschelman } 781c24bd37SHong Zhang 7958cf0668SJed Brown static PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_LLCondensed(Mat A,Mat B,PetscReal fill,Mat *C) 80b561aa9dSHong Zhang { 81b561aa9dSHong Zhang PetscErrorCode ierr; 82b561aa9dSHong Zhang Mat_SeqAIJ *a =(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c; 83971236abSHong Zhang PetscInt *ai=a->i,*bi=b->i,*ci,*cj; 84c1ba5575SJed Brown PetscInt am =A->rmap->N,bn=B->cmap->N,bm=B->rmap->N; 85b561aa9dSHong Zhang PetscReal afill; 86eca6b86aSHong Zhang PetscInt i,j,anzi,brow,bnzj,cnzi,*bj,*aj,*lnk,ndouble=0,Crmax; 87eca6b86aSHong Zhang PetscTable ta; 88fb908031SHong Zhang PetscBT lnkbt; 890298fd71SBarry Smith PetscFreeSpaceList free_space=NULL,current_space=NULL; 90b561aa9dSHong Zhang 91b561aa9dSHong Zhang PetscFunctionBegin; 92bd958071SHong Zhang /* Get ci and cj */ 93bd958071SHong Zhang /*---------------*/ 94fb908031SHong Zhang /* Allocate ci array, arrays for fill computation and */ 95fb908031SHong Zhang /* free space for accumulating nonzero column info */ 96854ce69bSBarry Smith ierr = PetscMalloc1(am+2,&ci);CHKERRQ(ierr); 97fb908031SHong Zhang ci[0] = 0; 98fb908031SHong Zhang 99fb908031SHong Zhang /* create and initialize a linked list */ 100c373ccc6SHong Zhang ierr = PetscTableCreate(bn,bn,&ta);CHKERRQ(ierr); 101c373ccc6SHong Zhang MatRowMergeMax_SeqAIJ(b,bm,ta); 102eca6b86aSHong Zhang ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr); 103eca6b86aSHong Zhang ierr = PetscTableDestroy(&ta);CHKERRQ(ierr); 104eca6b86aSHong Zhang 105eca6b86aSHong Zhang ierr = PetscLLCondensedCreate(Crmax,bn,&lnk,&lnkbt);CHKERRQ(ierr); 106fb908031SHong Zhang 107fb908031SHong Zhang /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */ 108f91af8c7SBarry Smith ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ai[am],bi[bm])),&free_space);CHKERRQ(ierr); 1092205254eSKarl Rupp 110fb908031SHong Zhang current_space = free_space; 111fb908031SHong Zhang 112fb908031SHong Zhang /* Determine ci and cj */ 113fb908031SHong Zhang for (i=0; i<am; i++) { 114fb908031SHong Zhang anzi = ai[i+1] - ai[i]; 115fb908031SHong Zhang aj = a->j + ai[i]; 116fb908031SHong Zhang for (j=0; j<anzi; j++) { 117fb908031SHong Zhang brow = aj[j]; 118fb908031SHong Zhang bnzj = bi[brow+1] - bi[brow]; 119fb908031SHong Zhang bj = b->j + bi[brow]; 120fb908031SHong Zhang /* add non-zero cols of B into the sorted linked list lnk */ 121fb908031SHong Zhang ierr = PetscLLCondensedAddSorted(bnzj,bj,lnk,lnkbt);CHKERRQ(ierr); 122fb908031SHong Zhang } 123fb908031SHong Zhang cnzi = lnk[0]; 124fb908031SHong Zhang 125fb908031SHong Zhang /* If free space is not available, make more free space */ 126fb908031SHong Zhang /* Double the amount of total space in the list */ 127fb908031SHong Zhang if (current_space->local_remaining<cnzi) { 128f91af8c7SBarry Smith ierr = PetscFreeSpaceGet(PetscIntSumTruncate(cnzi,current_space->total_array_size),¤t_space);CHKERRQ(ierr); 129fb908031SHong Zhang ndouble++; 130fb908031SHong Zhang } 131fb908031SHong Zhang 132fb908031SHong Zhang /* Copy data into free space, then initialize lnk */ 133fb908031SHong Zhang ierr = PetscLLCondensedClean(bn,cnzi,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 1342205254eSKarl Rupp 135fb908031SHong Zhang current_space->array += cnzi; 136fb908031SHong Zhang current_space->local_used += cnzi; 137fb908031SHong Zhang current_space->local_remaining -= cnzi; 1382205254eSKarl Rupp 139fb908031SHong Zhang ci[i+1] = ci[i] + cnzi; 140fb908031SHong Zhang } 141fb908031SHong Zhang 142fb908031SHong Zhang /* Column indices are in the list of free space */ 143fb908031SHong Zhang /* Allocate space for cj, initialize cj, and */ 144fb908031SHong Zhang /* destroy list of free space and other temporary array(s) */ 145854ce69bSBarry Smith ierr = PetscMalloc1(ci[am]+1,&cj);CHKERRQ(ierr); 146fb908031SHong Zhang ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr); 147fb908031SHong Zhang ierr = PetscLLCondensedDestroy(lnk,lnkbt);CHKERRQ(ierr); 148b561aa9dSHong Zhang 14926be0446SHong Zhang /* put together the new symbolic matrix */ 150ce94432eSBarry Smith ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,NULL,C);CHKERRQ(ierr); 15133d57670SJed Brown ierr = MatSetBlockSizesFromMats(*C,A,B);CHKERRQ(ierr); 15258c24d83SHong Zhang 15358c24d83SHong Zhang /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 15458c24d83SHong Zhang /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */ 15558c24d83SHong Zhang c = (Mat_SeqAIJ*)((*C)->data); 156aa1aec99SHong Zhang c->free_a = PETSC_FALSE; 157e6b907acSBarry Smith c->free_ij = PETSC_TRUE; 15858c24d83SHong Zhang c->nonew = 0; 159002e8affSHong Zhang (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ; /* fast, needs non-scalable O(bn) array 'abdense' */ 1600b7e3e3dSHong Zhang 1617212da7cSHong Zhang /* set MatInfo */ 1627212da7cSHong Zhang afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5; 163f2b054eeSHong Zhang if (afill < 1.0) afill = 1.0; 1647212da7cSHong Zhang c->maxnz = ci[am]; 1657212da7cSHong Zhang c->nz = ci[am]; 166fb908031SHong Zhang (*C)->info.mallocs = ndouble; 1677212da7cSHong Zhang (*C)->info.fill_ratio_given = fill; 1687212da7cSHong Zhang (*C)->info.fill_ratio_needed = afill; 1697212da7cSHong Zhang 1707212da7cSHong Zhang #if defined(PETSC_USE_INFO) 1717212da7cSHong Zhang if (ci[am]) { 17257622a8eSBarry Smith ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %g needed %g.\n",ndouble,(double)fill,(double)afill);CHKERRQ(ierr); 17357622a8eSBarry Smith ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr); 174f2b054eeSHong Zhang } else { 175f2b054eeSHong Zhang ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr); 176be0fcf8dSHong Zhang } 177f2b054eeSHong Zhang #endif 17858c24d83SHong Zhang PetscFunctionReturn(0); 17958c24d83SHong Zhang } 180d50806bdSBarry Smith 181dfbe8321SBarry Smith PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C) 182d50806bdSBarry Smith { 183dfbe8321SBarry Smith PetscErrorCode ierr; 184d13dce4bSSatish Balay PetscLogDouble flops=0.0; 185d50806bdSBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 186d50806bdSBarry Smith Mat_SeqAIJ *b = (Mat_SeqAIJ*)B->data; 187d50806bdSBarry Smith Mat_SeqAIJ *c = (Mat_SeqAIJ*)C->data; 18838baddfdSBarry Smith PetscInt *ai =a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j; 189c58ca2e3SHong Zhang PetscInt am =A->rmap->n,cm=C->rmap->n; 190519eb980SHong Zhang PetscInt i,j,k,anzi,bnzi,cnzi,brow; 191aa1aec99SHong Zhang PetscScalar *aa=a->a,*ba=b->a,*baj,*ca,valtmp; 192aa1aec99SHong Zhang PetscScalar *ab_dense; 193d50806bdSBarry Smith 194d50806bdSBarry Smith PetscFunctionBegin; 195aa1aec99SHong Zhang if (!c->a) { /* first call of MatMatMultNumeric_SeqAIJ_SeqAIJ, allocate ca and matmult_abdense */ 196854ce69bSBarry Smith ierr = PetscMalloc1(ci[cm]+1,&ca);CHKERRQ(ierr); 197aa1aec99SHong Zhang c->a = ca; 198aa1aec99SHong Zhang c->free_a = PETSC_TRUE; 199aa1aec99SHong Zhang } else { 200aa1aec99SHong Zhang ca = c->a; 201d90d86d0SMatthew G. Knepley } 202d90d86d0SMatthew G. Knepley if (!c->matmult_abdense) { 2031795a4d1SJed Brown ierr = PetscCalloc1(B->cmap->N,&ab_dense);CHKERRQ(ierr); 204d90d86d0SMatthew G. Knepley c->matmult_abdense = ab_dense; 205d90d86d0SMatthew G. Knepley } else { 206aa1aec99SHong Zhang ab_dense = c->matmult_abdense; 207aa1aec99SHong Zhang } 208aa1aec99SHong Zhang 209c124e916SHong Zhang /* clean old values in C */ 210c124e916SHong Zhang ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr); 211d50806bdSBarry Smith /* Traverse A row-wise. */ 212d50806bdSBarry Smith /* Build the ith row in C by summing over nonzero columns in A, */ 213d50806bdSBarry Smith /* the rows of B corresponding to nonzeros of A. */ 214d50806bdSBarry Smith for (i=0; i<am; i++) { 215d50806bdSBarry Smith anzi = ai[i+1] - ai[i]; 216d50806bdSBarry Smith for (j=0; j<anzi; j++) { 217519eb980SHong Zhang brow = aj[j]; 218d50806bdSBarry Smith bnzi = bi[brow+1] - bi[brow]; 219d50806bdSBarry Smith bjj = bj + bi[brow]; 220d50806bdSBarry Smith baj = ba + bi[brow]; 221519eb980SHong Zhang /* perform dense axpy */ 22236ec6d2dSHong Zhang valtmp = aa[j]; 223519eb980SHong Zhang for (k=0; k<bnzi; k++) { 22436ec6d2dSHong Zhang ab_dense[bjj[k]] += valtmp*baj[k]; 225519eb980SHong Zhang } 226519eb980SHong Zhang flops += 2*bnzi; 227519eb980SHong Zhang } 228c58ca2e3SHong Zhang aj += anzi; aa += anzi; 229c58ca2e3SHong Zhang 230c58ca2e3SHong Zhang cnzi = ci[i+1] - ci[i]; 231519eb980SHong Zhang for (k=0; k<cnzi; k++) { 232519eb980SHong Zhang ca[k] += ab_dense[cj[k]]; 233519eb980SHong Zhang ab_dense[cj[k]] = 0.0; /* zero ab_dense */ 234519eb980SHong Zhang } 235519eb980SHong Zhang flops += cnzi; 236519eb980SHong Zhang cj += cnzi; ca += cnzi; 237519eb980SHong Zhang } 238c58ca2e3SHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 239c58ca2e3SHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 240c58ca2e3SHong Zhang ierr = PetscLogFlops(flops);CHKERRQ(ierr); 241c58ca2e3SHong Zhang PetscFunctionReturn(0); 242c58ca2e3SHong Zhang } 243c58ca2e3SHong Zhang 24425023028SHong Zhang PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(Mat A,Mat B,Mat C) 245c58ca2e3SHong Zhang { 246c58ca2e3SHong Zhang PetscErrorCode ierr; 247c58ca2e3SHong Zhang PetscLogDouble flops=0.0; 248c58ca2e3SHong Zhang Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 249c58ca2e3SHong Zhang Mat_SeqAIJ *b = (Mat_SeqAIJ*)B->data; 250c58ca2e3SHong Zhang Mat_SeqAIJ *c = (Mat_SeqAIJ*)C->data; 251c58ca2e3SHong Zhang PetscInt *ai = a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j; 252c58ca2e3SHong Zhang PetscInt am = A->rmap->N,cm=C->rmap->N; 253c58ca2e3SHong Zhang PetscInt i,j,k,anzi,bnzi,cnzi,brow; 25436ec6d2dSHong Zhang PetscScalar *aa=a->a,*ba=b->a,*baj,*ca=c->a,valtmp; 255c58ca2e3SHong Zhang PetscInt nextb; 256c58ca2e3SHong Zhang 257c58ca2e3SHong Zhang PetscFunctionBegin; 258cf742fccSHong Zhang if (!ca) { /* first call of MatMatMultNumeric_SeqAIJ_SeqAIJ, allocate ca and matmult_abdense */ 259cf742fccSHong Zhang ierr = PetscMalloc1(ci[cm]+1,&ca);CHKERRQ(ierr); 260cf742fccSHong Zhang c->a = ca; 261cf742fccSHong Zhang c->free_a = PETSC_TRUE; 262cf742fccSHong Zhang } 263cf742fccSHong Zhang 264c58ca2e3SHong Zhang /* clean old values in C */ 265c58ca2e3SHong Zhang ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr); 266c58ca2e3SHong Zhang /* Traverse A row-wise. */ 267c58ca2e3SHong Zhang /* Build the ith row in C by summing over nonzero columns in A, */ 268c58ca2e3SHong Zhang /* the rows of B corresponding to nonzeros of A. */ 269519eb980SHong Zhang for (i=0; i<am; i++) { 270519eb980SHong Zhang anzi = ai[i+1] - ai[i]; 271519eb980SHong Zhang cnzi = ci[i+1] - ci[i]; 272519eb980SHong Zhang for (j=0; j<anzi; j++) { 273519eb980SHong Zhang brow = aj[j]; 274519eb980SHong Zhang bnzi = bi[brow+1] - bi[brow]; 275519eb980SHong Zhang bjj = bj + bi[brow]; 276519eb980SHong Zhang baj = ba + bi[brow]; 277519eb980SHong Zhang /* perform sparse axpy */ 27836ec6d2dSHong Zhang valtmp = aa[j]; 279c124e916SHong Zhang nextb = 0; 280c124e916SHong Zhang for (k=0; nextb<bnzi; k++) { 281c124e916SHong Zhang if (cj[k] == bjj[nextb]) { /* ccol == bcol */ 28236ec6d2dSHong Zhang ca[k] += valtmp*baj[nextb++]; 283c124e916SHong Zhang } 284d50806bdSBarry Smith } 285d50806bdSBarry Smith flops += 2*bnzi; 286d50806bdSBarry Smith } 287519eb980SHong Zhang aj += anzi; aa += anzi; 288519eb980SHong Zhang cj += cnzi; ca += cnzi; 289d50806bdSBarry Smith } 290c58ca2e3SHong Zhang 291716bacf3SKris Buschelman ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 292716bacf3SKris Buschelman ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 293d50806bdSBarry Smith ierr = PetscLogFlops(flops);CHKERRQ(ierr); 294d50806bdSBarry Smith PetscFunctionReturn(0); 295d50806bdSBarry Smith } 296bc011b1eSHong Zhang 2973c50cad2SHong Zhang PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast(Mat A,Mat B,PetscReal fill,Mat *C) 29825296bd5SBarry Smith { 29925296bd5SBarry Smith PetscErrorCode ierr; 30025296bd5SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c; 30125296bd5SBarry Smith PetscInt *ai = a->i,*bi=b->i,*ci,*cj; 3023c50cad2SHong Zhang PetscInt am = A->rmap->N,bn=B->cmap->N,bm=B->rmap->N; 30325296bd5SBarry Smith MatScalar *ca; 30425296bd5SBarry Smith PetscReal afill; 305eca6b86aSHong Zhang PetscInt i,j,anzi,brow,bnzj,cnzi,*bj,*aj,*lnk,ndouble=0,Crmax; 306eca6b86aSHong Zhang PetscTable ta; 3070298fd71SBarry Smith PetscFreeSpaceList free_space=NULL,current_space=NULL; 30825296bd5SBarry Smith 30925296bd5SBarry Smith PetscFunctionBegin; 3103c50cad2SHong Zhang /* Get ci and cj - same as MatMatMultSymbolic_SeqAIJ_SeqAIJ except using PetscLLxxx_fast() */ 3113c50cad2SHong Zhang /*-----------------------------------------------------------------------------------------*/ 3123c50cad2SHong Zhang /* Allocate arrays for fill computation and free space for accumulating nonzero column */ 313854ce69bSBarry Smith ierr = PetscMalloc1(am+2,&ci);CHKERRQ(ierr); 3143c50cad2SHong Zhang ci[0] = 0; 3153c50cad2SHong Zhang 3163c50cad2SHong Zhang /* create and initialize a linked list */ 317c373ccc6SHong Zhang ierr = PetscTableCreate(bn,bn,&ta);CHKERRQ(ierr); 318c373ccc6SHong Zhang MatRowMergeMax_SeqAIJ(b,bm,ta); 319eca6b86aSHong Zhang ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr); 320eca6b86aSHong Zhang ierr = PetscTableDestroy(&ta);CHKERRQ(ierr); 321eca6b86aSHong Zhang 322eca6b86aSHong Zhang ierr = PetscLLCondensedCreate_fast(Crmax,&lnk);CHKERRQ(ierr); 3233c50cad2SHong Zhang 3243c50cad2SHong Zhang /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */ 325f91af8c7SBarry Smith ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ai[am],bi[bm])),&free_space);CHKERRQ(ierr); 3263c50cad2SHong Zhang current_space = free_space; 3273c50cad2SHong Zhang 3283c50cad2SHong Zhang /* Determine ci and cj */ 3293c50cad2SHong Zhang for (i=0; i<am; i++) { 3303c50cad2SHong Zhang anzi = ai[i+1] - ai[i]; 3313c50cad2SHong Zhang aj = a->j + ai[i]; 3323c50cad2SHong Zhang for (j=0; j<anzi; j++) { 3333c50cad2SHong Zhang brow = aj[j]; 3343c50cad2SHong Zhang bnzj = bi[brow+1] - bi[brow]; 3353c50cad2SHong Zhang bj = b->j + bi[brow]; 3363c50cad2SHong Zhang /* add non-zero cols of B into the sorted linked list lnk */ 3373c50cad2SHong Zhang ierr = PetscLLCondensedAddSorted_fast(bnzj,bj,lnk);CHKERRQ(ierr); 3383c50cad2SHong Zhang } 3393c50cad2SHong Zhang cnzi = lnk[1]; 3403c50cad2SHong Zhang 3413c50cad2SHong Zhang /* If free space is not available, make more free space */ 3423c50cad2SHong Zhang /* Double the amount of total space in the list */ 3433c50cad2SHong Zhang if (current_space->local_remaining<cnzi) { 344f91af8c7SBarry Smith ierr = PetscFreeSpaceGet(PetscIntSumTruncate(cnzi,current_space->total_array_size),¤t_space);CHKERRQ(ierr); 3453c50cad2SHong Zhang ndouble++; 3463c50cad2SHong Zhang } 3473c50cad2SHong Zhang 3483c50cad2SHong Zhang /* Copy data into free space, then initialize lnk */ 3493c50cad2SHong Zhang ierr = PetscLLCondensedClean_fast(cnzi,current_space->array,lnk);CHKERRQ(ierr); 3502205254eSKarl Rupp 3513c50cad2SHong Zhang current_space->array += cnzi; 3523c50cad2SHong Zhang current_space->local_used += cnzi; 3533c50cad2SHong Zhang current_space->local_remaining -= cnzi; 3542205254eSKarl Rupp 3553c50cad2SHong Zhang ci[i+1] = ci[i] + cnzi; 3563c50cad2SHong Zhang } 3573c50cad2SHong Zhang 3583c50cad2SHong Zhang /* Column indices are in the list of free space */ 3593c50cad2SHong Zhang /* Allocate space for cj, initialize cj, and */ 3603c50cad2SHong Zhang /* destroy list of free space and other temporary array(s) */ 361854ce69bSBarry Smith ierr = PetscMalloc1(ci[am]+1,&cj);CHKERRQ(ierr); 3623c50cad2SHong Zhang ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr); 3633c50cad2SHong Zhang ierr = PetscLLCondensedDestroy_fast(lnk);CHKERRQ(ierr); 36425296bd5SBarry Smith 36525296bd5SBarry Smith /* Allocate space for ca */ 366854ce69bSBarry Smith ierr = PetscMalloc1(ci[am]+1,&ca);CHKERRQ(ierr); 36725296bd5SBarry Smith ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr); 36825296bd5SBarry Smith 36925296bd5SBarry Smith /* put together the new symbolic matrix */ 370ce94432eSBarry Smith ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,ca,C);CHKERRQ(ierr); 37133d57670SJed Brown ierr = MatSetBlockSizesFromMats(*C,A,B);CHKERRQ(ierr); 37225296bd5SBarry Smith 37325296bd5SBarry Smith /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 37425296bd5SBarry Smith /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */ 37525296bd5SBarry Smith c = (Mat_SeqAIJ*)((*C)->data); 37625296bd5SBarry Smith c->free_a = PETSC_TRUE; 37725296bd5SBarry Smith c->free_ij = PETSC_TRUE; 37825296bd5SBarry Smith c->nonew = 0; 3792205254eSKarl Rupp 38025296bd5SBarry Smith (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable; /* slower, less memory */ 38125296bd5SBarry Smith 38225296bd5SBarry Smith /* set MatInfo */ 38325296bd5SBarry Smith afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5; 38425296bd5SBarry Smith if (afill < 1.0) afill = 1.0; 38525296bd5SBarry Smith c->maxnz = ci[am]; 38625296bd5SBarry Smith c->nz = ci[am]; 3873c50cad2SHong Zhang (*C)->info.mallocs = ndouble; 38825296bd5SBarry Smith (*C)->info.fill_ratio_given = fill; 38925296bd5SBarry Smith (*C)->info.fill_ratio_needed = afill; 39025296bd5SBarry Smith 39125296bd5SBarry Smith #if defined(PETSC_USE_INFO) 39225296bd5SBarry Smith if (ci[am]) { 39357622a8eSBarry Smith ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %g needed %g.\n",ndouble,(double)fill,(double)afill);CHKERRQ(ierr); 39457622a8eSBarry Smith ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr); 39525296bd5SBarry Smith } else { 39625296bd5SBarry Smith ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr); 39725296bd5SBarry Smith } 39825296bd5SBarry Smith #endif 39925296bd5SBarry Smith PetscFunctionReturn(0); 40025296bd5SBarry Smith } 40125296bd5SBarry Smith 40225296bd5SBarry Smith 40325023028SHong Zhang PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(Mat A,Mat B,PetscReal fill,Mat *C) 404e9e4536cSHong Zhang { 405e9e4536cSHong Zhang PetscErrorCode ierr; 406e9e4536cSHong Zhang Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c; 407bf9555e6SHong Zhang PetscInt *ai = a->i,*bi=b->i,*ci,*cj; 40825c41797SHong Zhang PetscInt am = A->rmap->N,bn=B->cmap->N,bm=B->rmap->N; 409e9e4536cSHong Zhang MatScalar *ca; 410e9e4536cSHong Zhang PetscReal afill; 411eca6b86aSHong Zhang PetscInt i,j,anzi,brow,bnzj,cnzi,*bj,*aj,*lnk,ndouble=0,Crmax; 412eca6b86aSHong Zhang PetscTable ta; 4130298fd71SBarry Smith PetscFreeSpaceList free_space=NULL,current_space=NULL; 414e9e4536cSHong Zhang 415e9e4536cSHong Zhang PetscFunctionBegin; 416bd958071SHong Zhang /* Get ci and cj - same as MatMatMultSymbolic_SeqAIJ_SeqAIJ except using PetscLLxxx_Scalalbe() */ 417bd958071SHong Zhang /*---------------------------------------------------------------------------------------------*/ 418bd958071SHong Zhang /* Allocate arrays for fill computation and free space for accumulating nonzero column */ 419854ce69bSBarry Smith ierr = PetscMalloc1(am+2,&ci);CHKERRQ(ierr); 420bd958071SHong Zhang ci[0] = 0; 421bd958071SHong Zhang 422bd958071SHong Zhang /* create and initialize a linked list */ 423c373ccc6SHong Zhang ierr = PetscTableCreate(bn,bn,&ta);CHKERRQ(ierr); 424c373ccc6SHong Zhang MatRowMergeMax_SeqAIJ(b,bm,ta); 425eca6b86aSHong Zhang ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr); 426eca6b86aSHong Zhang ierr = PetscTableDestroy(&ta);CHKERRQ(ierr); 427eca6b86aSHong Zhang ierr = PetscLLCondensedCreate_Scalable(Crmax,&lnk);CHKERRQ(ierr); 428bd958071SHong Zhang 429bd958071SHong Zhang /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */ 430f91af8c7SBarry Smith ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ai[am],bi[bm])),&free_space);CHKERRQ(ierr); 431bd958071SHong Zhang current_space = free_space; 432bd958071SHong Zhang 433bd958071SHong Zhang /* Determine ci and cj */ 434bd958071SHong Zhang for (i=0; i<am; i++) { 435bd958071SHong Zhang anzi = ai[i+1] - ai[i]; 436bd958071SHong Zhang aj = a->j + ai[i]; 437bd958071SHong Zhang for (j=0; j<anzi; j++) { 438bd958071SHong Zhang brow = aj[j]; 439bd958071SHong Zhang bnzj = bi[brow+1] - bi[brow]; 440bd958071SHong Zhang bj = b->j + bi[brow]; 441bd958071SHong Zhang /* add non-zero cols of B into the sorted linked list lnk */ 442bd958071SHong Zhang ierr = PetscLLCondensedAddSorted_Scalable(bnzj,bj,lnk);CHKERRQ(ierr); 443bd958071SHong Zhang } 444bd958071SHong Zhang cnzi = lnk[0]; 445bd958071SHong Zhang 446bd958071SHong Zhang /* If free space is not available, make more free space */ 447bd958071SHong Zhang /* Double the amount of total space in the list */ 448bd958071SHong Zhang if (current_space->local_remaining<cnzi) { 449f91af8c7SBarry Smith ierr = PetscFreeSpaceGet(PetscIntSumTruncate(cnzi,current_space->total_array_size),¤t_space);CHKERRQ(ierr); 450bd958071SHong Zhang ndouble++; 451bd958071SHong Zhang } 452bd958071SHong Zhang 453bd958071SHong Zhang /* Copy data into free space, then initialize lnk */ 454bd958071SHong Zhang ierr = PetscLLCondensedClean_Scalable(cnzi,current_space->array,lnk);CHKERRQ(ierr); 4552205254eSKarl Rupp 456bd958071SHong Zhang current_space->array += cnzi; 457bd958071SHong Zhang current_space->local_used += cnzi; 458bd958071SHong Zhang current_space->local_remaining -= cnzi; 4592205254eSKarl Rupp 460bd958071SHong Zhang ci[i+1] = ci[i] + cnzi; 461bd958071SHong Zhang } 462bd958071SHong Zhang 463bd958071SHong Zhang /* Column indices are in the list of free space */ 464bd958071SHong Zhang /* Allocate space for cj, initialize cj, and */ 465bd958071SHong Zhang /* destroy list of free space and other temporary array(s) */ 466854ce69bSBarry Smith ierr = PetscMalloc1(ci[am]+1,&cj);CHKERRQ(ierr); 467bd958071SHong Zhang ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr); 468bd958071SHong Zhang ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr); 469e9e4536cSHong Zhang 470e9e4536cSHong Zhang /* Allocate space for ca */ 471bd958071SHong Zhang /*-----------------------*/ 472854ce69bSBarry Smith ierr = PetscMalloc1(ci[am]+1,&ca);CHKERRQ(ierr); 473e9e4536cSHong Zhang ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr); 474e9e4536cSHong Zhang 475e9e4536cSHong Zhang /* put together the new symbolic matrix */ 476ce94432eSBarry Smith ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,ca,C);CHKERRQ(ierr); 47733d57670SJed Brown ierr = MatSetBlockSizesFromMats(*C,A,B);CHKERRQ(ierr); 478e9e4536cSHong Zhang 479e9e4536cSHong Zhang /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 480e9e4536cSHong Zhang /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */ 481e9e4536cSHong Zhang c = (Mat_SeqAIJ*)((*C)->data); 482e9e4536cSHong Zhang c->free_a = PETSC_TRUE; 483e9e4536cSHong Zhang c->free_ij = PETSC_TRUE; 484e9e4536cSHong Zhang c->nonew = 0; 4852205254eSKarl Rupp 48625023028SHong Zhang (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable; /* slower, less memory */ 487e9e4536cSHong Zhang 488e9e4536cSHong Zhang /* set MatInfo */ 489e9e4536cSHong Zhang afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5; 490e9e4536cSHong Zhang if (afill < 1.0) afill = 1.0; 491e9e4536cSHong Zhang c->maxnz = ci[am]; 492e9e4536cSHong Zhang c->nz = ci[am]; 493bd958071SHong Zhang (*C)->info.mallocs = ndouble; 494e9e4536cSHong Zhang (*C)->info.fill_ratio_given = fill; 495e9e4536cSHong Zhang (*C)->info.fill_ratio_needed = afill; 496e9e4536cSHong Zhang 497e9e4536cSHong Zhang #if defined(PETSC_USE_INFO) 498e9e4536cSHong Zhang if (ci[am]) { 49957622a8eSBarry Smith ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %g needed %g.\n",ndouble,(double)fill,(double)afill);CHKERRQ(ierr); 50057622a8eSBarry Smith ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr); 501e9e4536cSHong Zhang } else { 502e9e4536cSHong Zhang ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr); 503e9e4536cSHong Zhang } 504e9e4536cSHong Zhang #endif 505e9e4536cSHong Zhang PetscFunctionReturn(0); 506e9e4536cSHong Zhang } 507e9e4536cSHong Zhang 5080ced3a2bSJed Brown PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Heap(Mat A,Mat B,PetscReal fill,Mat *C) 5090ced3a2bSJed Brown { 5100ced3a2bSJed Brown PetscErrorCode ierr; 5110ced3a2bSJed Brown Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c; 5120ced3a2bSJed Brown const PetscInt *ai=a->i,*bi=b->i,*aj=a->j,*bj=b->j; 5130ced3a2bSJed Brown PetscInt *ci,*cj,*bb; 5140ced3a2bSJed Brown PetscInt am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N; 5150ced3a2bSJed Brown PetscReal afill; 5160ced3a2bSJed Brown PetscInt i,j,col,ndouble = 0; 5170298fd71SBarry Smith PetscFreeSpaceList free_space=NULL,current_space=NULL; 5180ced3a2bSJed Brown PetscHeap h; 5190ced3a2bSJed Brown 5200ced3a2bSJed Brown PetscFunctionBegin; 521cd093f1dSJed Brown /* Get ci and cj - by merging sorted rows using a heap */ 5220ced3a2bSJed Brown /*---------------------------------------------------------------------------------------------*/ 5230ced3a2bSJed Brown /* Allocate arrays for fill computation and free space for accumulating nonzero column */ 524854ce69bSBarry Smith ierr = PetscMalloc1(am+2,&ci);CHKERRQ(ierr); 5250ced3a2bSJed Brown ci[0] = 0; 5260ced3a2bSJed Brown 5270ced3a2bSJed Brown /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */ 528f91af8c7SBarry Smith ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ai[am],bi[bm])),&free_space);CHKERRQ(ierr); 5290ced3a2bSJed Brown current_space = free_space; 5300ced3a2bSJed Brown 5310ced3a2bSJed Brown ierr = PetscHeapCreate(a->rmax,&h);CHKERRQ(ierr); 532785e854fSJed Brown ierr = PetscMalloc1(a->rmax,&bb);CHKERRQ(ierr); 5330ced3a2bSJed Brown 5340ced3a2bSJed Brown /* Determine ci and cj */ 5350ced3a2bSJed Brown for (i=0; i<am; i++) { 5360ced3a2bSJed 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 */ 5370ced3a2bSJed Brown const PetscInt *acol = aj + ai[i]; /* column indices of nonzero entries in this row */ 5380ced3a2bSJed Brown ci[i+1] = ci[i]; 5390ced3a2bSJed Brown /* Populate the min heap */ 5400ced3a2bSJed Brown for (j=0; j<anzi; j++) { 5410ced3a2bSJed Brown bb[j] = bi[acol[j]]; /* bb points at the start of the row */ 5420ced3a2bSJed Brown if (bb[j] < bi[acol[j]+1]) { /* Add if row is nonempty */ 5430ced3a2bSJed Brown ierr = PetscHeapAdd(h,j,bj[bb[j]++]);CHKERRQ(ierr); 5440ced3a2bSJed Brown } 5450ced3a2bSJed Brown } 5460ced3a2bSJed Brown /* Pick off the min element, adding it to free space */ 5470ced3a2bSJed Brown ierr = PetscHeapPop(h,&j,&col);CHKERRQ(ierr); 5480ced3a2bSJed Brown while (j >= 0) { 5490ced3a2bSJed Brown if (current_space->local_remaining < 1) { /* double the size, but don't exceed 16 MiB */ 550f91af8c7SBarry Smith ierr = PetscFreeSpaceGet(PetscMin(PetscIntMultTruncate(2,current_space->total_array_size),16 << 20),¤t_space);CHKERRQ(ierr); 5510ced3a2bSJed Brown ndouble++; 5520ced3a2bSJed Brown } 5530ced3a2bSJed Brown *(current_space->array++) = col; 5540ced3a2bSJed Brown current_space->local_used++; 5550ced3a2bSJed Brown current_space->local_remaining--; 5560ced3a2bSJed Brown ci[i+1]++; 5570ced3a2bSJed Brown 5580ced3a2bSJed Brown /* stash if anything else remains in this row of B */ 5590ced3a2bSJed Brown if (bb[j] < bi[acol[j]+1]) {ierr = PetscHeapStash(h,j,bj[bb[j]++]);CHKERRQ(ierr);} 5600ced3a2bSJed Brown while (1) { /* pop and stash any other rows of B that also had an entry in this column */ 5610ced3a2bSJed Brown PetscInt j2,col2; 5620ced3a2bSJed Brown ierr = PetscHeapPeek(h,&j2,&col2);CHKERRQ(ierr); 5630ced3a2bSJed Brown if (col2 != col) break; 5640ced3a2bSJed Brown ierr = PetscHeapPop(h,&j2,&col2);CHKERRQ(ierr); 5650ced3a2bSJed Brown if (bb[j2] < bi[acol[j2]+1]) {ierr = PetscHeapStash(h,j2,bj[bb[j2]++]);CHKERRQ(ierr);} 5660ced3a2bSJed Brown } 5670ced3a2bSJed Brown /* Put any stashed elements back into the min heap */ 5680ced3a2bSJed Brown ierr = PetscHeapUnstash(h);CHKERRQ(ierr); 5690ced3a2bSJed Brown ierr = PetscHeapPop(h,&j,&col);CHKERRQ(ierr); 5700ced3a2bSJed Brown } 5710ced3a2bSJed Brown } 5720ced3a2bSJed Brown ierr = PetscFree(bb);CHKERRQ(ierr); 5730ced3a2bSJed Brown ierr = PetscHeapDestroy(&h);CHKERRQ(ierr); 5740ced3a2bSJed Brown 5750ced3a2bSJed Brown /* Column indices are in the list of free space */ 5760ced3a2bSJed Brown /* Allocate space for cj, initialize cj, and */ 5770ced3a2bSJed Brown /* destroy list of free space and other temporary array(s) */ 578785e854fSJed Brown ierr = PetscMalloc1(ci[am],&cj);CHKERRQ(ierr); 5790ced3a2bSJed Brown ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr); 5800ced3a2bSJed Brown 5810ced3a2bSJed Brown /* put together the new symbolic matrix */ 582ce94432eSBarry Smith ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,NULL,C);CHKERRQ(ierr); 58333d57670SJed Brown ierr = MatSetBlockSizesFromMats(*C,A,B);CHKERRQ(ierr); 5840ced3a2bSJed Brown 5850ced3a2bSJed Brown /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 5860ced3a2bSJed Brown /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */ 5870ced3a2bSJed Brown c = (Mat_SeqAIJ*)((*C)->data); 5880ced3a2bSJed Brown c->free_a = PETSC_TRUE; 5890ced3a2bSJed Brown c->free_ij = PETSC_TRUE; 5900ced3a2bSJed Brown c->nonew = 0; 59126fbe8dcSKarl Rupp 59289d95c1aSJed Brown (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ; 5930ced3a2bSJed Brown 5940ced3a2bSJed Brown /* set MatInfo */ 5950ced3a2bSJed Brown afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5; 5960ced3a2bSJed Brown if (afill < 1.0) afill = 1.0; 5970ced3a2bSJed Brown c->maxnz = ci[am]; 5980ced3a2bSJed Brown c->nz = ci[am]; 5990ced3a2bSJed Brown (*C)->info.mallocs = ndouble; 6000ced3a2bSJed Brown (*C)->info.fill_ratio_given = fill; 6010ced3a2bSJed Brown (*C)->info.fill_ratio_needed = afill; 6020ced3a2bSJed Brown 6030ced3a2bSJed Brown #if defined(PETSC_USE_INFO) 6040ced3a2bSJed Brown if (ci[am]) { 60557622a8eSBarry Smith ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %g needed %g.\n",ndouble,(double)fill,(double)afill);CHKERRQ(ierr); 60657622a8eSBarry Smith ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr); 6070ced3a2bSJed Brown } else { 6080ced3a2bSJed Brown ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr); 6090ced3a2bSJed Brown } 6100ced3a2bSJed Brown #endif 6110ced3a2bSJed Brown PetscFunctionReturn(0); 6120ced3a2bSJed Brown } 613e9e4536cSHong Zhang 6148a07c6f1SJed Brown PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_BTHeap(Mat A,Mat B,PetscReal fill,Mat *C) 6158a07c6f1SJed Brown { 6168a07c6f1SJed Brown PetscErrorCode ierr; 6178a07c6f1SJed Brown Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c; 6188a07c6f1SJed Brown const PetscInt *ai = a->i,*bi=b->i,*aj=a->j,*bj=b->j; 6198a07c6f1SJed Brown PetscInt *ci,*cj,*bb; 6208a07c6f1SJed Brown PetscInt am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N; 6218a07c6f1SJed Brown PetscReal afill; 6228a07c6f1SJed Brown PetscInt i,j,col,ndouble = 0; 6230298fd71SBarry Smith PetscFreeSpaceList free_space=NULL,current_space=NULL; 6248a07c6f1SJed Brown PetscHeap h; 6258a07c6f1SJed Brown PetscBT bt; 6268a07c6f1SJed Brown 6278a07c6f1SJed Brown PetscFunctionBegin; 628cd093f1dSJed Brown /* Get ci and cj - using a heap for the sorted rows, but use BT so that each index is only added once */ 6298a07c6f1SJed Brown /*---------------------------------------------------------------------------------------------*/ 6308a07c6f1SJed Brown /* Allocate arrays for fill computation and free space for accumulating nonzero column */ 631854ce69bSBarry Smith ierr = PetscMalloc1(am+2,&ci);CHKERRQ(ierr); 6328a07c6f1SJed Brown ci[0] = 0; 6338a07c6f1SJed Brown 6348a07c6f1SJed Brown /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */ 635f91af8c7SBarry Smith ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ai[am],bi[bm])),&free_space);CHKERRQ(ierr); 6362205254eSKarl Rupp 6378a07c6f1SJed Brown current_space = free_space; 6388a07c6f1SJed Brown 6398a07c6f1SJed Brown ierr = PetscHeapCreate(a->rmax,&h);CHKERRQ(ierr); 640785e854fSJed Brown ierr = PetscMalloc1(a->rmax,&bb);CHKERRQ(ierr); 6418a07c6f1SJed Brown ierr = PetscBTCreate(bn,&bt);CHKERRQ(ierr); 6428a07c6f1SJed Brown 6438a07c6f1SJed Brown /* Determine ci and cj */ 6448a07c6f1SJed Brown for (i=0; i<am; i++) { 6458a07c6f1SJed 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 */ 6468a07c6f1SJed Brown const PetscInt *acol = aj + ai[i]; /* column indices of nonzero entries in this row */ 6478a07c6f1SJed Brown const PetscInt *fptr = current_space->array; /* Save beginning of the row so we can clear the BT later */ 6488a07c6f1SJed Brown ci[i+1] = ci[i]; 6498a07c6f1SJed Brown /* Populate the min heap */ 6508a07c6f1SJed Brown for (j=0; j<anzi; j++) { 6518a07c6f1SJed Brown PetscInt brow = acol[j]; 6528a07c6f1SJed Brown for (bb[j] = bi[brow]; bb[j] < bi[brow+1]; bb[j]++) { 6538a07c6f1SJed Brown PetscInt bcol = bj[bb[j]]; 6548a07c6f1SJed Brown if (!PetscBTLookupSet(bt,bcol)) { /* new entry */ 6558a07c6f1SJed Brown ierr = PetscHeapAdd(h,j,bcol);CHKERRQ(ierr); 6568a07c6f1SJed Brown bb[j]++; 6578a07c6f1SJed Brown break; 6588a07c6f1SJed Brown } 6598a07c6f1SJed Brown } 6608a07c6f1SJed Brown } 6618a07c6f1SJed Brown /* Pick off the min element, adding it to free space */ 6628a07c6f1SJed Brown ierr = PetscHeapPop(h,&j,&col);CHKERRQ(ierr); 6638a07c6f1SJed Brown while (j >= 0) { 6648a07c6f1SJed Brown if (current_space->local_remaining < 1) { /* double the size, but don't exceed 16 MiB */ 6650298fd71SBarry Smith fptr = NULL; /* need PetscBTMemzero */ 666f91af8c7SBarry Smith ierr = PetscFreeSpaceGet(PetscMin(PetscIntMultTruncate(2,current_space->total_array_size),16 << 20),¤t_space);CHKERRQ(ierr); 6678a07c6f1SJed Brown ndouble++; 6688a07c6f1SJed Brown } 6698a07c6f1SJed Brown *(current_space->array++) = col; 6708a07c6f1SJed Brown current_space->local_used++; 6718a07c6f1SJed Brown current_space->local_remaining--; 6728a07c6f1SJed Brown ci[i+1]++; 6738a07c6f1SJed Brown 6748a07c6f1SJed Brown /* stash if anything else remains in this row of B */ 6758a07c6f1SJed Brown for (; bb[j] < bi[acol[j]+1]; bb[j]++) { 6768a07c6f1SJed Brown PetscInt bcol = bj[bb[j]]; 6778a07c6f1SJed Brown if (!PetscBTLookupSet(bt,bcol)) { /* new entry */ 6788a07c6f1SJed Brown ierr = PetscHeapAdd(h,j,bcol);CHKERRQ(ierr); 6798a07c6f1SJed Brown bb[j]++; 6808a07c6f1SJed Brown break; 6818a07c6f1SJed Brown } 6828a07c6f1SJed Brown } 6838a07c6f1SJed Brown ierr = PetscHeapPop(h,&j,&col);CHKERRQ(ierr); 6848a07c6f1SJed Brown } 6858a07c6f1SJed Brown if (fptr) { /* Clear the bits for this row */ 6868a07c6f1SJed Brown for (; fptr<current_space->array; fptr++) {ierr = PetscBTClear(bt,*fptr);CHKERRQ(ierr);} 6878a07c6f1SJed Brown } else { /* We reallocated so we don't remember (easily) how to clear only the bits we changed */ 6888a07c6f1SJed Brown ierr = PetscBTMemzero(bn,bt);CHKERRQ(ierr); 6898a07c6f1SJed Brown } 6908a07c6f1SJed Brown } 6918a07c6f1SJed Brown ierr = PetscFree(bb);CHKERRQ(ierr); 6928a07c6f1SJed Brown ierr = PetscHeapDestroy(&h);CHKERRQ(ierr); 6938a07c6f1SJed Brown ierr = PetscBTDestroy(&bt);CHKERRQ(ierr); 6948a07c6f1SJed Brown 6958a07c6f1SJed Brown /* Column indices are in the list of free space */ 6968a07c6f1SJed Brown /* Allocate space for cj, initialize cj, and */ 6978a07c6f1SJed Brown /* destroy list of free space and other temporary array(s) */ 698785e854fSJed Brown ierr = PetscMalloc1(ci[am],&cj);CHKERRQ(ierr); 6998a07c6f1SJed Brown ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr); 7008a07c6f1SJed Brown 7018a07c6f1SJed Brown /* put together the new symbolic matrix */ 702ce94432eSBarry Smith ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,NULL,C);CHKERRQ(ierr); 70333d57670SJed Brown ierr = MatSetBlockSizesFromMats(*C,A,B);CHKERRQ(ierr); 7048a07c6f1SJed Brown 7058a07c6f1SJed Brown /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 7068a07c6f1SJed Brown /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */ 7078a07c6f1SJed Brown c = (Mat_SeqAIJ*)((*C)->data); 7088a07c6f1SJed Brown c->free_a = PETSC_TRUE; 7098a07c6f1SJed Brown c->free_ij = PETSC_TRUE; 7108a07c6f1SJed Brown c->nonew = 0; 71126fbe8dcSKarl Rupp 71289d95c1aSJed Brown (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ; 7138a07c6f1SJed Brown 7148a07c6f1SJed Brown /* set MatInfo */ 7158a07c6f1SJed Brown afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5; 7168a07c6f1SJed Brown if (afill < 1.0) afill = 1.0; 7178a07c6f1SJed Brown c->maxnz = ci[am]; 7188a07c6f1SJed Brown c->nz = ci[am]; 7198a07c6f1SJed Brown (*C)->info.mallocs = ndouble; 7208a07c6f1SJed Brown (*C)->info.fill_ratio_given = fill; 7218a07c6f1SJed Brown (*C)->info.fill_ratio_needed = afill; 7228a07c6f1SJed Brown 7238a07c6f1SJed Brown #if defined(PETSC_USE_INFO) 7248a07c6f1SJed Brown if (ci[am]) { 72557622a8eSBarry Smith ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %g needed %g.\n",ndouble,(double)fill,(double)afill);CHKERRQ(ierr); 72657622a8eSBarry Smith ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr); 7278a07c6f1SJed Brown } else { 7288a07c6f1SJed Brown ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr); 7298a07c6f1SJed Brown } 7308a07c6f1SJed Brown #endif 7318a07c6f1SJed Brown PetscFunctionReturn(0); 7328a07c6f1SJed Brown } 7338a07c6f1SJed Brown 734cd093f1dSJed Brown /* concatenate unique entries and then sort */ 73558cf0668SJed Brown PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C) 736cd093f1dSJed Brown { 737cd093f1dSJed Brown PetscErrorCode ierr; 738cd093f1dSJed Brown Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c; 739cd093f1dSJed Brown const PetscInt *ai = a->i,*bi=b->i,*aj=a->j,*bj=b->j; 740cd093f1dSJed Brown PetscInt *ci,*cj; 741cd093f1dSJed Brown PetscInt am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N; 742cd093f1dSJed Brown PetscReal afill; 743cd093f1dSJed Brown PetscInt i,j,ndouble = 0; 744cd093f1dSJed Brown PetscSegBuffer seg,segrow; 745cd093f1dSJed Brown char *seen; 746cd093f1dSJed Brown 747cd093f1dSJed Brown PetscFunctionBegin; 748854ce69bSBarry Smith ierr = PetscMalloc1(am+1,&ci);CHKERRQ(ierr); 749cd093f1dSJed Brown ci[0] = 0; 750cd093f1dSJed Brown 751cd093f1dSJed Brown /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */ 752cd093f1dSJed Brown ierr = PetscSegBufferCreate(sizeof(PetscInt),(PetscInt)(fill*(ai[am]+bi[bm])),&seg);CHKERRQ(ierr); 753cd093f1dSJed Brown ierr = PetscSegBufferCreate(sizeof(PetscInt),100,&segrow);CHKERRQ(ierr); 754785e854fSJed Brown ierr = PetscMalloc1(bn,&seen);CHKERRQ(ierr); 755cd093f1dSJed Brown ierr = PetscMemzero(seen,bn*sizeof(char));CHKERRQ(ierr); 756cd093f1dSJed Brown 757cd093f1dSJed Brown /* Determine ci and cj */ 758cd093f1dSJed Brown for (i=0; i<am; i++) { 759cd093f1dSJed 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 */ 760cd093f1dSJed Brown const PetscInt *acol = aj + ai[i]; /* column indices of nonzero entries in this row */ 761cd093f1dSJed Brown PetscInt packlen = 0,*PETSC_RESTRICT crow; 762cd093f1dSJed Brown /* Pack segrow */ 763cd093f1dSJed Brown for (j=0; j<anzi; j++) { 764cd093f1dSJed Brown PetscInt brow = acol[j],bjstart = bi[brow],bjend = bi[brow+1],k; 765cd093f1dSJed Brown for (k=bjstart; k<bjend; k++) { 766cd093f1dSJed Brown PetscInt bcol = bj[k]; 767cd093f1dSJed Brown if (!seen[bcol]) { /* new entry */ 768cd093f1dSJed Brown PetscInt *PETSC_RESTRICT slot; 769cd093f1dSJed Brown ierr = PetscSegBufferGetInts(segrow,1,&slot);CHKERRQ(ierr); 770cd093f1dSJed Brown *slot = bcol; 771cd093f1dSJed Brown seen[bcol] = 1; 772cd093f1dSJed Brown packlen++; 773cd093f1dSJed Brown } 774cd093f1dSJed Brown } 775cd093f1dSJed Brown } 776cd093f1dSJed Brown ierr = PetscSegBufferGetInts(seg,packlen,&crow);CHKERRQ(ierr); 777cd093f1dSJed Brown ierr = PetscSegBufferExtractTo(segrow,crow);CHKERRQ(ierr); 778cd093f1dSJed Brown ierr = PetscSortInt(packlen,crow);CHKERRQ(ierr); 779cd093f1dSJed Brown ci[i+1] = ci[i] + packlen; 780cd093f1dSJed Brown for (j=0; j<packlen; j++) seen[crow[j]] = 0; 781cd093f1dSJed Brown } 782cd093f1dSJed Brown ierr = PetscSegBufferDestroy(&segrow);CHKERRQ(ierr); 783cd093f1dSJed Brown ierr = PetscFree(seen);CHKERRQ(ierr); 784cd093f1dSJed Brown 785cd093f1dSJed Brown /* Column indices are in the segmented buffer */ 786cd093f1dSJed Brown ierr = PetscSegBufferExtractAlloc(seg,&cj);CHKERRQ(ierr); 787cd093f1dSJed Brown ierr = PetscSegBufferDestroy(&seg);CHKERRQ(ierr); 788cd093f1dSJed Brown 789cd093f1dSJed Brown /* put together the new symbolic matrix */ 790cd093f1dSJed Brown ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,NULL,C);CHKERRQ(ierr); 79133d57670SJed Brown ierr = MatSetBlockSizesFromMats(*C,A,B);CHKERRQ(ierr); 792cd093f1dSJed Brown 793cd093f1dSJed Brown /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 794cd093f1dSJed Brown /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */ 795cd093f1dSJed Brown c = (Mat_SeqAIJ*)((*C)->data); 796cd093f1dSJed Brown c->free_a = PETSC_TRUE; 797cd093f1dSJed Brown c->free_ij = PETSC_TRUE; 798cd093f1dSJed Brown c->nonew = 0; 799cd093f1dSJed Brown 800cd093f1dSJed Brown (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ; 801cd093f1dSJed Brown 802cd093f1dSJed Brown /* set MatInfo */ 803cd093f1dSJed Brown afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5; 804cd093f1dSJed Brown if (afill < 1.0) afill = 1.0; 805cd093f1dSJed Brown c->maxnz = ci[am]; 806cd093f1dSJed Brown c->nz = ci[am]; 807cd093f1dSJed Brown (*C)->info.mallocs = ndouble; 808cd093f1dSJed Brown (*C)->info.fill_ratio_given = fill; 809cd093f1dSJed Brown (*C)->info.fill_ratio_needed = afill; 810cd093f1dSJed Brown 811cd093f1dSJed Brown #if defined(PETSC_USE_INFO) 812cd093f1dSJed Brown if (ci[am]) { 81357622a8eSBarry Smith ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %g needed %g.\n",ndouble,(double)fill,(double)afill);CHKERRQ(ierr); 81457622a8eSBarry Smith ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr); 815cd093f1dSJed Brown } else { 816cd093f1dSJed Brown ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr); 817cd093f1dSJed Brown } 818cd093f1dSJed Brown #endif 819cd093f1dSJed Brown PetscFunctionReturn(0); 820cd093f1dSJed Brown } 821cd093f1dSJed Brown 822d2853540SHong Zhang /* This routine is not used. Should be removed! */ 8236fc122caSHong Zhang PetscErrorCode MatMatTransposeMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 8245df89d91SHong Zhang { 825bc011b1eSHong Zhang PetscErrorCode ierr; 826bc011b1eSHong Zhang 827bc011b1eSHong Zhang PetscFunctionBegin; 828bc011b1eSHong Zhang if (scall == MAT_INITIAL_MATRIX) { 8293ff4c91cSHong Zhang ierr = PetscLogEventBegin(MAT_MatTransposeMultSymbolic,A,B,0,0);CHKERRQ(ierr); 8306fc122caSHong Zhang ierr = MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr); 8313ff4c91cSHong Zhang ierr = PetscLogEventEnd(MAT_MatTransposeMultSymbolic,A,B,0,0);CHKERRQ(ierr); 832bc011b1eSHong Zhang } 8333ff4c91cSHong Zhang ierr = PetscLogEventBegin(MAT_MatTransposeMultNumeric,A,B,0,0);CHKERRQ(ierr); 8346fc122caSHong Zhang ierr = MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr); 8353ff4c91cSHong Zhang ierr = PetscLogEventEnd(MAT_MatTransposeMultNumeric,A,B,0,0);CHKERRQ(ierr); 836bc011b1eSHong Zhang PetscFunctionReturn(0); 837bc011b1eSHong Zhang } 838bc011b1eSHong Zhang 8392128a86cSHong Zhang PetscErrorCode MatDestroy_SeqAIJ_MatMatMultTrans(Mat A) 8402128a86cSHong Zhang { 8412128a86cSHong Zhang PetscErrorCode ierr; 8424c7df5ccSHong Zhang Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data; 84340192850SHong Zhang Mat_MatMatTransMult *abt=a->abt; 8442128a86cSHong Zhang 8452128a86cSHong Zhang PetscFunctionBegin; 84640192850SHong Zhang ierr = (abt->destroy)(A);CHKERRQ(ierr); 84740192850SHong Zhang ierr = MatTransposeColoringDestroy(&abt->matcoloring);CHKERRQ(ierr); 84840192850SHong Zhang ierr = MatDestroy(&abt->Bt_den);CHKERRQ(ierr); 84940192850SHong Zhang ierr = MatDestroy(&abt->ABt_den);CHKERRQ(ierr); 85040192850SHong Zhang ierr = PetscFree(abt);CHKERRQ(ierr); 8512128a86cSHong Zhang PetscFunctionReturn(0); 8522128a86cSHong Zhang } 8532128a86cSHong Zhang 8546fc122caSHong Zhang PetscErrorCode MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C) 855bc011b1eSHong Zhang { 8565df89d91SHong Zhang PetscErrorCode ierr; 85781d82fe4SHong Zhang Mat Bt; 85881d82fe4SHong Zhang PetscInt *bti,*btj; 85940192850SHong Zhang Mat_MatMatTransMult *abt; 8604c7df5ccSHong Zhang Mat_SeqAIJ *c; 861d2853540SHong Zhang 86281d82fe4SHong Zhang PetscFunctionBegin; 86381d82fe4SHong Zhang /* create symbolic Bt */ 86481d82fe4SHong Zhang ierr = MatGetSymbolicTranspose_SeqAIJ(B,&bti,&btj);CHKERRQ(ierr); 8650298fd71SBarry Smith ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,B->cmap->n,B->rmap->n,bti,btj,NULL,&Bt);CHKERRQ(ierr); 86633d57670SJed Brown ierr = MatSetBlockSizes(Bt,PetscAbs(A->cmap->bs),PetscAbs(B->cmap->bs));CHKERRQ(ierr); 86781d82fe4SHong Zhang 86881d82fe4SHong Zhang /* get symbolic C=A*Bt */ 86981d82fe4SHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,Bt,fill,C);CHKERRQ(ierr); 87081d82fe4SHong Zhang 8712128a86cSHong Zhang /* create a supporting struct for reuse intermidiate dense matrices with matcoloring */ 872b00a9115SJed Brown ierr = PetscNew(&abt);CHKERRQ(ierr); 8734c7df5ccSHong Zhang c = (Mat_SeqAIJ*)(*C)->data; 87440192850SHong Zhang c->abt = abt; 8752128a86cSHong Zhang 87640192850SHong Zhang abt->usecoloring = PETSC_FALSE; 87740192850SHong Zhang abt->destroy = (*C)->ops->destroy; 8782128a86cSHong Zhang (*C)->ops->destroy = MatDestroy_SeqAIJ_MatMatMultTrans; 8792128a86cSHong Zhang 880c5929fdfSBarry Smith ierr = PetscOptionsGetBool(((PetscObject)A)->options,NULL,"-matmattransmult_color",&abt->usecoloring,NULL);CHKERRQ(ierr); 88140192850SHong Zhang if (abt->usecoloring) { 882b9af6bddSHong Zhang /* Create MatTransposeColoring from symbolic C=A*B^T */ 883b9af6bddSHong Zhang MatTransposeColoring matcoloring; 884335efc43SPeter Brune MatColoring coloring; 8852128a86cSHong Zhang ISColoring iscoloring; 8862128a86cSHong Zhang Mat Bt_dense,C_dense; 8874d478ae7SHong Zhang Mat_SeqAIJ *c=(Mat_SeqAIJ*)(*C)->data; 8884d478ae7SHong Zhang /* inode causes memory problem, don't know why */ 8894d478ae7SHong Zhang if (c->inode.use) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MAT_USE_INODES is not supported. Use '-mat_no_inode'"); 890e8354b3eSHong Zhang 891335efc43SPeter Brune ierr = MatColoringCreate(*C,&coloring);CHKERRQ(ierr); 892335efc43SPeter Brune ierr = MatColoringSetDistance(coloring,2);CHKERRQ(ierr); 893335efc43SPeter Brune ierr = MatColoringSetType(coloring,MATCOLORINGSL);CHKERRQ(ierr); 894335efc43SPeter Brune ierr = MatColoringSetFromOptions(coloring);CHKERRQ(ierr); 895335efc43SPeter Brune ierr = MatColoringApply(coloring,&iscoloring);CHKERRQ(ierr); 896335efc43SPeter Brune ierr = MatColoringDestroy(&coloring);CHKERRQ(ierr); 897b9af6bddSHong Zhang ierr = MatTransposeColoringCreate(*C,iscoloring,&matcoloring);CHKERRQ(ierr); 8982205254eSKarl Rupp 89940192850SHong Zhang abt->matcoloring = matcoloring; 9002205254eSKarl Rupp 9012128a86cSHong Zhang ierr = ISColoringDestroy(&iscoloring);CHKERRQ(ierr); 9022128a86cSHong Zhang 9032128a86cSHong Zhang /* Create Bt_dense and C_dense = A*Bt_dense */ 9042128a86cSHong Zhang ierr = MatCreate(PETSC_COMM_SELF,&Bt_dense);CHKERRQ(ierr); 9052128a86cSHong Zhang ierr = MatSetSizes(Bt_dense,A->cmap->n,matcoloring->ncolors,A->cmap->n,matcoloring->ncolors);CHKERRQ(ierr); 9062128a86cSHong Zhang ierr = MatSetType(Bt_dense,MATSEQDENSE);CHKERRQ(ierr); 9070298fd71SBarry Smith ierr = MatSeqDenseSetPreallocation(Bt_dense,NULL);CHKERRQ(ierr); 9082205254eSKarl Rupp 9092128a86cSHong Zhang Bt_dense->assembled = PETSC_TRUE; 91040192850SHong Zhang abt->Bt_den = Bt_dense; 9112128a86cSHong Zhang 9122128a86cSHong Zhang ierr = MatCreate(PETSC_COMM_SELF,&C_dense);CHKERRQ(ierr); 9132128a86cSHong Zhang ierr = MatSetSizes(C_dense,A->rmap->n,matcoloring->ncolors,A->rmap->n,matcoloring->ncolors);CHKERRQ(ierr); 9142128a86cSHong Zhang ierr = MatSetType(C_dense,MATSEQDENSE);CHKERRQ(ierr); 9150298fd71SBarry Smith ierr = MatSeqDenseSetPreallocation(C_dense,NULL);CHKERRQ(ierr); 9162205254eSKarl Rupp 9172128a86cSHong Zhang Bt_dense->assembled = PETSC_TRUE; 91840192850SHong Zhang abt->ABt_den = C_dense; 919f94ccd6cSHong Zhang 920f94ccd6cSHong Zhang #if defined(PETSC_USE_INFO) 9211ce71dffSSatish Balay { 922f94ccd6cSHong Zhang Mat_SeqAIJ *c = (Mat_SeqAIJ*)(*C)->data; 923c40ebe3bSHong Zhang ierr = PetscInfo7(*C,"Use coloring of C=A*B^T; B^T: %D %D, Bt_dense: %D,%D; Cnz %D / (cm*ncolors %D) = %g\n",B->cmap->n,B->rmap->n,Bt_dense->rmap->n,Bt_dense->cmap->n,c->nz,A->rmap->n*matcoloring->ncolors,(PetscReal)(c->nz)/(A->rmap->n*matcoloring->ncolors));CHKERRQ(ierr); 9241ce71dffSSatish Balay } 925f94ccd6cSHong Zhang #endif 9262128a86cSHong Zhang } 927e99be685SHong Zhang /* clean up */ 928e99be685SHong Zhang ierr = MatDestroy(&Bt);CHKERRQ(ierr); 929e99be685SHong Zhang ierr = MatRestoreSymbolicTranspose_SeqAIJ(B,&bti,&btj);CHKERRQ(ierr); 9305df89d91SHong Zhang PetscFunctionReturn(0); 9315df89d91SHong Zhang } 9325df89d91SHong Zhang 9336fc122caSHong Zhang PetscErrorCode MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C) 9345df89d91SHong Zhang { 9355df89d91SHong Zhang PetscErrorCode ierr; 9365df89d91SHong Zhang Mat_SeqAIJ *a =(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data; 937e2cac8caSJed Brown PetscInt *ai =a->i,*aj=a->j,*bi=b->i,*bj=b->j,anzi,bnzj,nexta,nextb,*acol,*bcol,brow; 93881d82fe4SHong Zhang PetscInt cm =C->rmap->n,*ci=c->i,*cj=c->j,i,j,cnzi,*ccol; 9395df89d91SHong Zhang PetscLogDouble flops=0.0; 940aa1aec99SHong Zhang MatScalar *aa =a->a,*aval,*ba=b->a,*bval,*ca,*cval; 94140192850SHong Zhang Mat_MatMatTransMult *abt = c->abt; 9425df89d91SHong Zhang 9435df89d91SHong Zhang PetscFunctionBegin; 94458ed3ceaSHong Zhang /* clear old values in C */ 94558ed3ceaSHong Zhang if (!c->a) { 946854ce69bSBarry Smith ierr = PetscMalloc1(ci[cm]+1,&ca);CHKERRQ(ierr); 94758ed3ceaSHong Zhang c->a = ca; 94858ed3ceaSHong Zhang c->free_a = PETSC_TRUE; 94958ed3ceaSHong Zhang } else { 95058ed3ceaSHong Zhang ca = c->a; 95158ed3ceaSHong Zhang } 95258ed3ceaSHong Zhang ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr); 95358ed3ceaSHong Zhang 95440192850SHong Zhang if (abt->usecoloring) { 95540192850SHong Zhang MatTransposeColoring matcoloring = abt->matcoloring; 95640192850SHong Zhang Mat Bt_dense,C_dense = abt->ABt_den; 957c8db22f6SHong Zhang 958b9af6bddSHong Zhang /* Get Bt_dense by Apply MatTransposeColoring to B */ 95940192850SHong Zhang Bt_dense = abt->Bt_den; 960b9af6bddSHong Zhang ierr = MatTransColoringApplySpToDen(matcoloring,B,Bt_dense);CHKERRQ(ierr); 961c8db22f6SHong Zhang 962c8db22f6SHong Zhang /* C_dense = A*Bt_dense */ 9632128a86cSHong Zhang ierr = MatMatMultNumeric_SeqAIJ_SeqDense(A,Bt_dense,C_dense);CHKERRQ(ierr); 964c8db22f6SHong Zhang 9652128a86cSHong Zhang /* Recover C from C_dense */ 966b9af6bddSHong Zhang ierr = MatTransColoringApplyDenToSp(matcoloring,C_dense,C);CHKERRQ(ierr); 967c8db22f6SHong Zhang PetscFunctionReturn(0); 968c8db22f6SHong Zhang } 969c8db22f6SHong Zhang 9701fa1209cSHong Zhang for (i=0; i<cm; i++) { 97181d82fe4SHong Zhang anzi = ai[i+1] - ai[i]; 97281d82fe4SHong Zhang acol = aj + ai[i]; 9736973769bSHong Zhang aval = aa + ai[i]; 9741fa1209cSHong Zhang cnzi = ci[i+1] - ci[i]; 9751fa1209cSHong Zhang ccol = cj + ci[i]; 9766973769bSHong Zhang cval = ca + ci[i]; 9771fa1209cSHong Zhang for (j=0; j<cnzi; j++) { 97881d82fe4SHong Zhang brow = ccol[j]; 97981d82fe4SHong Zhang bnzj = bi[brow+1] - bi[brow]; 98081d82fe4SHong Zhang bcol = bj + bi[brow]; 9816973769bSHong Zhang bval = ba + bi[brow]; 9826973769bSHong Zhang 98381d82fe4SHong Zhang /* perform sparse inner-product c(i,j)=A[i,:]*B[j,:]^T */ 98481d82fe4SHong Zhang nexta = 0; nextb = 0; 98581d82fe4SHong Zhang while (nexta<anzi && nextb<bnzj) { 9867b6d5e96SMark Adams while (nexta < anzi && acol[nexta] < bcol[nextb]) nexta++; 98781d82fe4SHong Zhang if (nexta == anzi) break; 9887b6d5e96SMark Adams while (nextb < bnzj && acol[nexta] > bcol[nextb]) nextb++; 98981d82fe4SHong Zhang if (nextb == bnzj) break; 99081d82fe4SHong Zhang if (acol[nexta] == bcol[nextb]) { 9916973769bSHong Zhang cval[j] += aval[nexta]*bval[nextb]; 99281d82fe4SHong Zhang nexta++; nextb++; 99381d82fe4SHong Zhang flops += 2; 9941fa1209cSHong Zhang } 9951fa1209cSHong Zhang } 99681d82fe4SHong Zhang } 99781d82fe4SHong Zhang } 99881d82fe4SHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 99981d82fe4SHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 100081d82fe4SHong Zhang ierr = PetscLogFlops(flops);CHKERRQ(ierr); 10015df89d91SHong Zhang PetscFunctionReturn(0); 10025df89d91SHong Zhang } 10035df89d91SHong Zhang 10046d373c3eSHong Zhang PetscErrorCode MatDestroy_SeqAIJ_MatTransMatMult(Mat A) 10056d373c3eSHong Zhang { 10066d373c3eSHong Zhang PetscErrorCode ierr; 10076d373c3eSHong Zhang Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 10086d373c3eSHong Zhang Mat_MatTransMatMult *atb = a->atb; 10096d373c3eSHong Zhang 10106d373c3eSHong Zhang PetscFunctionBegin; 10116d373c3eSHong Zhang ierr = MatDestroy(&atb->At);CHKERRQ(ierr); 10126d373c3eSHong Zhang ierr = (atb->destroy)(A);CHKERRQ(ierr); 10136d373c3eSHong Zhang ierr = PetscFree(atb);CHKERRQ(ierr); 10146d373c3eSHong Zhang PetscFunctionReturn(0); 10156d373c3eSHong Zhang } 10166d373c3eSHong Zhang 10170adebc6cSBarry Smith PetscErrorCode MatTransposeMatMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 10180adebc6cSBarry Smith { 10195df89d91SHong Zhang PetscErrorCode ierr; 10206d373c3eSHong Zhang const char *algTypes[2] = {"matmatmult","outerproduct"}; 10216d373c3eSHong Zhang PetscInt alg=0; /* set default algorithm */ 10226d373c3eSHong Zhang Mat At; 10236d373c3eSHong Zhang Mat_MatTransMatMult *atb; 10246d373c3eSHong Zhang Mat_SeqAIJ *c; 10255df89d91SHong Zhang 10265df89d91SHong Zhang PetscFunctionBegin; 10275df89d91SHong Zhang if (scall == MAT_INITIAL_MATRIX) { 10286d373c3eSHong Zhang ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr); 10296d373c3eSHong Zhang PetscOptionsObject->alreadyprinted = PETSC_FALSE; /* a hack to ensure the option shows in '-help' */ 10306d373c3eSHong Zhang ierr = PetscOptionsEList("-mattransposematmult_via","Algorithmic approach","MatTransposeMatMult",algTypes,2,algTypes[0],&alg,NULL);CHKERRQ(ierr); 10316d373c3eSHong Zhang ierr = PetscOptionsEnd();CHKERRQ(ierr); 10326d373c3eSHong Zhang 10336d373c3eSHong Zhang switch (alg) { 10346d373c3eSHong Zhang case 1: 103575648e8dSHong Zhang ierr = MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr); 10366d373c3eSHong Zhang break; 10376d373c3eSHong Zhang default: 10386d373c3eSHong Zhang ierr = PetscNew(&atb);CHKERRQ(ierr); 10396d373c3eSHong Zhang ierr = MatTranspose_SeqAIJ(A,MAT_INITIAL_MATRIX,&At);CHKERRQ(ierr); 10406d373c3eSHong Zhang ierr = MatMatMult_SeqAIJ_SeqAIJ(At,B,MAT_INITIAL_MATRIX,fill,C);CHKERRQ(ierr); 10416d373c3eSHong Zhang 1042618cf492SHong Zhang c = (Mat_SeqAIJ*)(*C)->data; 10436d373c3eSHong Zhang c->atb = atb; 10446d373c3eSHong Zhang atb->At = At; 10456d373c3eSHong Zhang atb->destroy = (*C)->ops->destroy; 10466d373c3eSHong Zhang (*C)->ops->destroy = MatDestroy_SeqAIJ_MatTransMatMult; 10476d373c3eSHong Zhang 10486d373c3eSHong Zhang break; 10495df89d91SHong Zhang } 10506d373c3eSHong Zhang } 10516d373c3eSHong Zhang if (alg) { 10526d373c3eSHong Zhang ierr = (*(*C)->ops->mattransposemultnumeric)(A,B,*C);CHKERRQ(ierr); 10536d373c3eSHong Zhang } else if (!alg && scall == MAT_REUSE_MATRIX) { 10546d373c3eSHong Zhang c = (Mat_SeqAIJ*)(*C)->data; 10556d373c3eSHong Zhang atb = c->atb; 10566d373c3eSHong Zhang At = atb->At; 10576d373c3eSHong Zhang ierr = MatTranspose_SeqAIJ(A,MAT_REUSE_MATRIX,&At);CHKERRQ(ierr); 10586d373c3eSHong Zhang ierr = MatMatMult_SeqAIJ_SeqAIJ(At,B,MAT_REUSE_MATRIX,fill,C);CHKERRQ(ierr); 10596d373c3eSHong Zhang } 10605df89d91SHong Zhang PetscFunctionReturn(0); 10615df89d91SHong Zhang } 10625df89d91SHong Zhang 106375648e8dSHong Zhang PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C) 10645df89d91SHong Zhang { 1065bc011b1eSHong Zhang PetscErrorCode ierr; 1066bc011b1eSHong Zhang Mat At; 106738baddfdSBarry Smith PetscInt *ati,*atj; 1068bc011b1eSHong Zhang 1069bc011b1eSHong Zhang PetscFunctionBegin; 1070bc011b1eSHong Zhang /* create symbolic At */ 1071bc011b1eSHong Zhang ierr = MatGetSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr); 10720298fd71SBarry Smith ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,A->cmap->n,A->rmap->n,ati,atj,NULL,&At);CHKERRQ(ierr); 107333d57670SJed Brown ierr = MatSetBlockSizes(At,PetscAbs(A->cmap->bs),PetscAbs(B->cmap->bs));CHKERRQ(ierr); 1074bc011b1eSHong Zhang 1075bc011b1eSHong Zhang /* get symbolic C=At*B */ 1076bc011b1eSHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(At,B,fill,C);CHKERRQ(ierr); 1077bc011b1eSHong Zhang 1078bc011b1eSHong Zhang /* clean up */ 10796bf464f9SBarry Smith ierr = MatDestroy(&At);CHKERRQ(ierr); 1080bc011b1eSHong Zhang ierr = MatRestoreSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr); 10816d373c3eSHong Zhang 10826d373c3eSHong Zhang (*C)->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ; 1083bc011b1eSHong Zhang PetscFunctionReturn(0); 1084bc011b1eSHong Zhang } 1085bc011b1eSHong Zhang 108675648e8dSHong Zhang PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C) 1087bc011b1eSHong Zhang { 1088bc011b1eSHong Zhang PetscErrorCode ierr; 10890fbc74f4SHong Zhang Mat_SeqAIJ *a =(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data; 1090d0f46423SBarry Smith PetscInt am =A->rmap->n,anzi,*ai=a->i,*aj=a->j,*bi=b->i,*bj,bnzi,nextb; 1091d0f46423SBarry Smith PetscInt cm =C->rmap->n,*ci=c->i,*cj=c->j,crow,*cjj,i,j,k; 1092d13dce4bSSatish Balay PetscLogDouble flops=0.0; 1093aa1aec99SHong Zhang MatScalar *aa =a->a,*ba,*ca,*caj; 1094bc011b1eSHong Zhang 1095bc011b1eSHong Zhang PetscFunctionBegin; 1096aa1aec99SHong Zhang if (!c->a) { 1097854ce69bSBarry Smith ierr = PetscMalloc1(ci[cm]+1,&ca);CHKERRQ(ierr); 10982205254eSKarl Rupp 1099aa1aec99SHong Zhang c->a = ca; 1100aa1aec99SHong Zhang c->free_a = PETSC_TRUE; 1101aa1aec99SHong Zhang } else { 1102aa1aec99SHong Zhang ca = c->a; 1103aa1aec99SHong Zhang } 1104bc011b1eSHong Zhang /* clear old values in C */ 1105bc011b1eSHong Zhang ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr); 1106bc011b1eSHong Zhang 1107bc011b1eSHong Zhang /* compute A^T*B using outer product (A^T)[:,i]*B[i,:] */ 1108bc011b1eSHong Zhang for (i=0; i<am; i++) { 1109bc011b1eSHong Zhang bj = b->j + bi[i]; 1110bc011b1eSHong Zhang ba = b->a + bi[i]; 1111bc011b1eSHong Zhang bnzi = bi[i+1] - bi[i]; 1112bc011b1eSHong Zhang anzi = ai[i+1] - ai[i]; 1113bc011b1eSHong Zhang for (j=0; j<anzi; j++) { 1114bc011b1eSHong Zhang nextb = 0; 11150fbc74f4SHong Zhang crow = *aj++; 1116bc011b1eSHong Zhang cjj = cj + ci[crow]; 1117bc011b1eSHong Zhang caj = ca + ci[crow]; 1118bc011b1eSHong Zhang /* perform sparse axpy operation. Note cjj includes bj. */ 1119bc011b1eSHong Zhang for (k=0; nextb<bnzi; k++) { 11200fbc74f4SHong Zhang if (cjj[k] == *(bj+nextb)) { /* ccol == bcol */ 11210fbc74f4SHong Zhang caj[k] += (*aa)*(*(ba+nextb)); 1122bc011b1eSHong Zhang nextb++; 1123bc011b1eSHong Zhang } 1124bc011b1eSHong Zhang } 1125bc011b1eSHong Zhang flops += 2*bnzi; 11260fbc74f4SHong Zhang aa++; 1127bc011b1eSHong Zhang } 1128bc011b1eSHong Zhang } 1129bc011b1eSHong Zhang 1130bc011b1eSHong Zhang /* Assemble the final matrix and clean up */ 1131bc011b1eSHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1132bc011b1eSHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1133bc011b1eSHong Zhang ierr = PetscLogFlops(flops);CHKERRQ(ierr); 1134bc011b1eSHong Zhang PetscFunctionReturn(0); 1135bc011b1eSHong Zhang } 11369123193aSHong Zhang 1137150d2497SBarry Smith PETSC_INTERN PetscErrorCode MatMatMult_SeqAIJ_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 11389123193aSHong Zhang { 11399123193aSHong Zhang PetscErrorCode ierr; 11409123193aSHong Zhang 11419123193aSHong Zhang PetscFunctionBegin; 11429123193aSHong Zhang if (scall == MAT_INITIAL_MATRIX) { 11433ff4c91cSHong Zhang ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr); 11449123193aSHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqDense(A,B,fill,C);CHKERRQ(ierr); 11453ff4c91cSHong Zhang ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr); 11469123193aSHong Zhang } 11473ff4c91cSHong Zhang ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr); 11489123193aSHong Zhang ierr = MatMatMultNumeric_SeqAIJ_SeqDense(A,B,*C);CHKERRQ(ierr); 11494614e006SHong Zhang ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr); 11509123193aSHong Zhang PetscFunctionReturn(0); 11519123193aSHong Zhang } 1152edd81eecSMatthew Knepley 11539123193aSHong Zhang PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C) 11549123193aSHong Zhang { 11559123193aSHong Zhang PetscErrorCode ierr; 11569123193aSHong Zhang 11579123193aSHong Zhang PetscFunctionBegin; 11585a586d82SBarry Smith ierr = MatMatMultSymbolic_SeqDense_SeqDense(A,B,0.0,C);CHKERRQ(ierr); 11592205254eSKarl Rupp 1160d73949e8SHong Zhang (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqDense; 11619123193aSHong Zhang PetscFunctionReturn(0); 11629123193aSHong Zhang } 11639123193aSHong Zhang 11649123193aSHong Zhang PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqDense(Mat A,Mat B,Mat C) 11659123193aSHong Zhang { 1166f32d5d43SBarry Smith Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data; 1167612438f1SStefano Zampini Mat_SeqDense *bd = (Mat_SeqDense*)B->data; 11689123193aSHong Zhang PetscErrorCode ierr; 1169daf57211SHong Zhang PetscScalar *c,*b,r1,r2,r3,r4,*c1,*c2,*c3,*c4,aatmp; 1170daf57211SHong Zhang const PetscScalar *aa,*b1,*b2,*b3,*b4; 1171daf57211SHong Zhang const PetscInt *aj; 1172612438f1SStefano Zampini PetscInt cm=C->rmap->n,cn=B->cmap->n,bm=bd->lda,am=A->rmap->n; 1173daf57211SHong Zhang PetscInt am4=4*am,bm4=4*bm,col,i,j,n,ajtmp; 11749123193aSHong Zhang 11759123193aSHong Zhang PetscFunctionBegin; 1176f32d5d43SBarry Smith if (!cm || !cn) PetscFunctionReturn(0); 1177612438f1SStefano Zampini if (B->rmap->n != 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,B->rmap->n); 1178e32f2f54SBarry 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); 1179e32f2f54SBarry 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); 1180612438f1SStefano Zampini b = bd->v; 11818c778c55SBarry Smith ierr = MatDenseGetArray(C,&c);CHKERRQ(ierr); 1182f32d5d43SBarry Smith b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm; 1183730858b9SHong Zhang c1 = c; c2 = c1 + am; c3 = c2 + am; c4 = c3 + am; 1184f32d5d43SBarry Smith for (col=0; col<cn-4; col += 4) { /* over columns of C */ 1185f32d5d43SBarry Smith for (i=0; i<am; i++) { /* over rows of C in those columns */ 1186f32d5d43SBarry Smith r1 = r2 = r3 = r4 = 0.0; 1187f32d5d43SBarry Smith n = a->i[i+1] - a->i[i]; 1188f32d5d43SBarry Smith aj = a->j + a->i[i]; 1189f32d5d43SBarry Smith aa = a->a + a->i[i]; 1190f32d5d43SBarry Smith for (j=0; j<n; j++) { 1191730858b9SHong Zhang aatmp = aa[j]; ajtmp = aj[j]; 1192730858b9SHong Zhang r1 += aatmp*b1[ajtmp]; 1193730858b9SHong Zhang r2 += aatmp*b2[ajtmp]; 1194730858b9SHong Zhang r3 += aatmp*b3[ajtmp]; 1195730858b9SHong Zhang r4 += aatmp*b4[ajtmp]; 11969123193aSHong Zhang } 1197730858b9SHong Zhang c1[i] = r1; 1198730858b9SHong Zhang c2[i] = r2; 1199730858b9SHong Zhang c3[i] = r3; 1200730858b9SHong Zhang c4[i] = r4; 1201f32d5d43SBarry Smith } 1202730858b9SHong Zhang b1 += bm4; b2 += bm4; b3 += bm4; b4 += bm4; 1203730858b9SHong Zhang c1 += am4; c2 += am4; c3 += am4; c4 += am4; 1204f32d5d43SBarry Smith } 1205f32d5d43SBarry Smith for (; col<cn; col++) { /* over extra columns of C */ 1206f32d5d43SBarry Smith for (i=0; i<am; i++) { /* over rows of C in those columns */ 1207f32d5d43SBarry Smith r1 = 0.0; 1208f32d5d43SBarry Smith n = a->i[i+1] - a->i[i]; 1209f32d5d43SBarry Smith aj = a->j + a->i[i]; 1210f32d5d43SBarry Smith aa = a->a + a->i[i]; 1211f32d5d43SBarry Smith for (j=0; j<n; j++) { 1212730858b9SHong Zhang r1 += aa[j]*b1[aj[j]]; 1213f32d5d43SBarry Smith } 1214730858b9SHong Zhang c1[i] = r1; 1215f32d5d43SBarry Smith } 1216f32d5d43SBarry Smith b1 += bm; 1217730858b9SHong Zhang c1 += am; 1218f32d5d43SBarry Smith } 1219dc0b31edSSatish Balay ierr = PetscLogFlops(cn*(2.0*a->nz));CHKERRQ(ierr); 12208c778c55SBarry Smith ierr = MatDenseRestoreArray(C,&c);CHKERRQ(ierr); 1221f32d5d43SBarry Smith ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1222f32d5d43SBarry Smith ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1223f32d5d43SBarry Smith PetscFunctionReturn(0); 1224f32d5d43SBarry Smith } 1225f32d5d43SBarry Smith 1226f32d5d43SBarry Smith /* 12274324174eSBarry Smith Note very similar to MatMult_SeqAIJ(), should generate both codes from same base 1228f32d5d43SBarry Smith */ 1229f32d5d43SBarry Smith PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat A,Mat B,Mat C) 1230f32d5d43SBarry Smith { 1231f32d5d43SBarry Smith Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data; 123290f5ea3eSStefano Zampini Mat_SeqDense *bd = (Mat_SeqDense*)B->data; 1233f32d5d43SBarry Smith PetscErrorCode ierr; 1234dd6ea824SBarry Smith PetscScalar *b,*c,r1,r2,r3,r4,*b1,*b2,*b3,*b4; 1235dd6ea824SBarry Smith MatScalar *aa; 123690f5ea3eSStefano Zampini PetscInt cm = C->rmap->n, cn=B->cmap->n, bm=bd->lda, col, i,j,n,*aj, am = A->rmap->n,*ii,arm; 12374324174eSBarry Smith PetscInt am2 = 2*am, am3 = 3*am, bm4 = 4*bm,colam,*ridx; 1238f32d5d43SBarry Smith 1239f32d5d43SBarry Smith PetscFunctionBegin; 1240f32d5d43SBarry Smith if (!cm || !cn) PetscFunctionReturn(0); 124190f5ea3eSStefano Zampini b = bd->v; 12428c778c55SBarry Smith ierr = MatDenseGetArray(C,&c);CHKERRQ(ierr); 1243f32d5d43SBarry Smith b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm; 12444324174eSBarry Smith 12454324174eSBarry Smith if (a->compressedrow.use) { /* use compressed row format */ 12464324174eSBarry Smith for (col=0; col<cn-4; col += 4) { /* over columns of C */ 12474324174eSBarry Smith colam = col*am; 12484324174eSBarry Smith arm = a->compressedrow.nrows; 12494324174eSBarry Smith ii = a->compressedrow.i; 12504324174eSBarry Smith ridx = a->compressedrow.rindex; 12514324174eSBarry Smith for (i=0; i<arm; i++) { /* over rows of C in those columns */ 12524324174eSBarry Smith r1 = r2 = r3 = r4 = 0.0; 12534324174eSBarry Smith n = ii[i+1] - ii[i]; 12544324174eSBarry Smith aj = a->j + ii[i]; 12554324174eSBarry Smith aa = a->a + ii[i]; 12564324174eSBarry Smith for (j=0; j<n; j++) { 12574324174eSBarry Smith r1 += (*aa)*b1[*aj]; 12584324174eSBarry Smith r2 += (*aa)*b2[*aj]; 12594324174eSBarry Smith r3 += (*aa)*b3[*aj]; 12604324174eSBarry Smith r4 += (*aa++)*b4[*aj++]; 12614324174eSBarry Smith } 12624324174eSBarry Smith c[colam + ridx[i]] += r1; 12634324174eSBarry Smith c[colam + am + ridx[i]] += r2; 12644324174eSBarry Smith c[colam + am2 + ridx[i]] += r3; 12654324174eSBarry Smith c[colam + am3 + ridx[i]] += r4; 12664324174eSBarry Smith } 12674324174eSBarry Smith b1 += bm4; 12684324174eSBarry Smith b2 += bm4; 12694324174eSBarry Smith b3 += bm4; 12704324174eSBarry Smith b4 += bm4; 12714324174eSBarry Smith } 12724324174eSBarry Smith for (; col<cn; col++) { /* over extra columns of C */ 12734324174eSBarry Smith colam = col*am; 12744324174eSBarry Smith arm = a->compressedrow.nrows; 12754324174eSBarry Smith ii = a->compressedrow.i; 12764324174eSBarry Smith ridx = a->compressedrow.rindex; 12774324174eSBarry Smith for (i=0; i<arm; i++) { /* over rows of C in those columns */ 12784324174eSBarry Smith r1 = 0.0; 12794324174eSBarry Smith n = ii[i+1] - ii[i]; 12804324174eSBarry Smith aj = a->j + ii[i]; 12814324174eSBarry Smith aa = a->a + ii[i]; 12824324174eSBarry Smith 12834324174eSBarry Smith for (j=0; j<n; j++) { 12844324174eSBarry Smith r1 += (*aa++)*b1[*aj++]; 12854324174eSBarry Smith } 1286a2ea699eSBarry Smith c[colam + ridx[i]] += r1; 12874324174eSBarry Smith } 12884324174eSBarry Smith b1 += bm; 12894324174eSBarry Smith } 12904324174eSBarry Smith } else { 1291f32d5d43SBarry Smith for (col=0; col<cn-4; col += 4) { /* over columns of C */ 1292f32d5d43SBarry Smith colam = col*am; 1293f32d5d43SBarry Smith for (i=0; i<am; i++) { /* over rows of C in those columns */ 1294f32d5d43SBarry Smith r1 = r2 = r3 = r4 = 0.0; 1295f32d5d43SBarry Smith n = a->i[i+1] - a->i[i]; 1296f32d5d43SBarry Smith aj = a->j + a->i[i]; 1297f32d5d43SBarry Smith aa = a->a + a->i[i]; 1298f32d5d43SBarry Smith for (j=0; j<n; j++) { 1299f32d5d43SBarry Smith r1 += (*aa)*b1[*aj]; 1300f32d5d43SBarry Smith r2 += (*aa)*b2[*aj]; 1301f32d5d43SBarry Smith r3 += (*aa)*b3[*aj]; 1302f32d5d43SBarry Smith r4 += (*aa++)*b4[*aj++]; 1303f32d5d43SBarry Smith } 1304f32d5d43SBarry Smith c[colam + i] += r1; 1305f32d5d43SBarry Smith c[colam + am + i] += r2; 1306f32d5d43SBarry Smith c[colam + am2 + i] += r3; 1307f32d5d43SBarry Smith c[colam + am3 + i] += r4; 1308f32d5d43SBarry Smith } 1309f32d5d43SBarry Smith b1 += bm4; 1310f32d5d43SBarry Smith b2 += bm4; 1311f32d5d43SBarry Smith b3 += bm4; 1312f32d5d43SBarry Smith b4 += bm4; 1313f32d5d43SBarry Smith } 1314f32d5d43SBarry Smith for (; col<cn; col++) { /* over extra columns of C */ 1315a2ea699eSBarry Smith colam = col*am; 1316f32d5d43SBarry Smith for (i=0; i<am; i++) { /* over rows of C in those columns */ 1317f32d5d43SBarry Smith r1 = 0.0; 1318f32d5d43SBarry Smith n = a->i[i+1] - a->i[i]; 1319f32d5d43SBarry Smith aj = a->j + a->i[i]; 1320f32d5d43SBarry Smith aa = a->a + a->i[i]; 1321f32d5d43SBarry Smith 1322f32d5d43SBarry Smith for (j=0; j<n; j++) { 1323f32d5d43SBarry Smith r1 += (*aa++)*b1[*aj++]; 1324f32d5d43SBarry Smith } 1325a2ea699eSBarry Smith c[colam + i] += r1; 1326f32d5d43SBarry Smith } 1327f32d5d43SBarry Smith b1 += bm; 1328f32d5d43SBarry Smith } 13294324174eSBarry Smith } 1330dc0b31edSSatish Balay ierr = PetscLogFlops(cn*2.0*a->nz);CHKERRQ(ierr); 13318c778c55SBarry Smith ierr = MatDenseRestoreArray(C,&c);CHKERRQ(ierr); 13329123193aSHong Zhang PetscFunctionReturn(0); 13339123193aSHong Zhang } 1334b1683b59SHong Zhang 1335b9af6bddSHong Zhang PetscErrorCode MatTransColoringApplySpToDen_SeqAIJ(MatTransposeColoring coloring,Mat B,Mat Btdense) 1336c8db22f6SHong Zhang { 1337c8db22f6SHong Zhang PetscErrorCode ierr; 13382f5f1f90SHong Zhang Mat_SeqAIJ *b = (Mat_SeqAIJ*)B->data; 13392f5f1f90SHong Zhang Mat_SeqDense *btdense = (Mat_SeqDense*)Btdense->data; 13402f5f1f90SHong Zhang PetscInt *bi = b->i,*bj=b->j; 13412f5f1f90SHong Zhang PetscInt m = Btdense->rmap->n,n=Btdense->cmap->n,j,k,l,col,anz,*btcol,brow,ncolumns; 13422f5f1f90SHong Zhang MatScalar *btval,*btval_den,*ba=b->a; 13432f5f1f90SHong Zhang PetscInt *columns=coloring->columns,*colorforcol=coloring->colorforcol,ncolors=coloring->ncolors; 1344c8db22f6SHong Zhang 1345c8db22f6SHong Zhang PetscFunctionBegin; 13462f5f1f90SHong Zhang btval_den=btdense->v; 13472f5f1f90SHong Zhang ierr = PetscMemzero(btval_den,(m*n)*sizeof(MatScalar));CHKERRQ(ierr); 13482f5f1f90SHong Zhang for (k=0; k<ncolors; k++) { 13492f5f1f90SHong Zhang ncolumns = coloring->ncolumns[k]; 13502f5f1f90SHong Zhang for (l=0; l<ncolumns; l++) { /* insert a row of B to a column of Btdense */ 1351d2853540SHong Zhang col = *(columns + colorforcol[k] + l); 13522f5f1f90SHong Zhang btcol = bj + bi[col]; 13532f5f1f90SHong Zhang btval = ba + bi[col]; 13542f5f1f90SHong Zhang anz = bi[col+1] - bi[col]; 1355d2853540SHong Zhang for (j=0; j<anz; j++) { 13562f5f1f90SHong Zhang brow = btcol[j]; 13572f5f1f90SHong Zhang btval_den[brow] = btval[j]; 1358c8db22f6SHong Zhang } 1359c8db22f6SHong Zhang } 13602f5f1f90SHong Zhang btval_den += m; 1361c8db22f6SHong Zhang } 1362c8db22f6SHong Zhang PetscFunctionReturn(0); 1363c8db22f6SHong Zhang } 1364c8db22f6SHong Zhang 1365b9af6bddSHong Zhang PetscErrorCode MatTransColoringApplyDenToSp_SeqAIJ(MatTransposeColoring matcoloring,Mat Cden,Mat Csp) 13668972f759SHong Zhang { 13678972f759SHong Zhang PetscErrorCode ierr; 1368b2d2b04fSHong Zhang Mat_SeqAIJ *csp = (Mat_SeqAIJ*)Csp->data; 1369077f23c2SHong Zhang PetscScalar *ca_den,*ca_den_ptr,*ca=csp->a; 1370f99a636bSHong Zhang PetscInt k,l,m=Cden->rmap->n,ncolors=matcoloring->ncolors; 1371e88777f2SHong Zhang PetscInt brows=matcoloring->brows,*den2sp=matcoloring->den2sp; 1372077f23c2SHong Zhang PetscInt nrows,*row,*idx; 1373077f23c2SHong Zhang PetscInt *rows=matcoloring->rows,*colorforrow=matcoloring->colorforrow; 13748972f759SHong Zhang 13758972f759SHong Zhang PetscFunctionBegin; 1376a3fe58edSHong Zhang ierr = MatDenseGetArray(Cden,&ca_den);CHKERRQ(ierr); 1377f99a636bSHong Zhang 1378077f23c2SHong Zhang if (brows > 0) { 1379077f23c2SHong Zhang PetscInt *lstart,row_end,row_start; 1380980ae229SHong Zhang lstart = matcoloring->lstart; 1381eeb4fd02SHong Zhang ierr = PetscMemzero(lstart,ncolors*sizeof(PetscInt));CHKERRQ(ierr); 1382eeb4fd02SHong Zhang 1383077f23c2SHong Zhang row_end = brows; 1384eeb4fd02SHong Zhang if (row_end > m) row_end = m; 1385077f23c2SHong Zhang for (row_start=0; row_start<m; row_start+=brows) { /* loop over row blocks of Csp */ 1386077f23c2SHong Zhang ca_den_ptr = ca_den; 1387980ae229SHong Zhang for (k=0; k<ncolors; k++) { /* loop over colors (columns of Cden) */ 1388eeb4fd02SHong Zhang nrows = matcoloring->nrows[k]; 1389eeb4fd02SHong Zhang row = rows + colorforrow[k]; 1390eeb4fd02SHong Zhang idx = den2sp + colorforrow[k]; 1391eeb4fd02SHong Zhang for (l=lstart[k]; l<nrows; l++) { 1392eeb4fd02SHong Zhang if (row[l] >= row_end) { 1393eeb4fd02SHong Zhang lstart[k] = l; 1394eeb4fd02SHong Zhang break; 1395eeb4fd02SHong Zhang } else { 1396077f23c2SHong Zhang ca[idx[l]] = ca_den_ptr[row[l]]; 1397eeb4fd02SHong Zhang } 1398eeb4fd02SHong Zhang } 1399077f23c2SHong Zhang ca_den_ptr += m; 1400eeb4fd02SHong Zhang } 1401077f23c2SHong Zhang row_end += brows; 1402eeb4fd02SHong Zhang if (row_end > m) row_end = m; 1403eeb4fd02SHong Zhang } 1404077f23c2SHong Zhang } else { /* non-blocked impl: loop over columns of Csp - slow if Csp is large */ 1405077f23c2SHong Zhang ca_den_ptr = ca_den; 1406077f23c2SHong Zhang for (k=0; k<ncolors; k++) { 1407077f23c2SHong Zhang nrows = matcoloring->nrows[k]; 1408077f23c2SHong Zhang row = rows + colorforrow[k]; 1409077f23c2SHong Zhang idx = den2sp + colorforrow[k]; 1410077f23c2SHong Zhang for (l=0; l<nrows; l++) { 1411077f23c2SHong Zhang ca[idx[l]] = ca_den_ptr[row[l]]; 1412077f23c2SHong Zhang } 1413077f23c2SHong Zhang ca_den_ptr += m; 1414077f23c2SHong Zhang } 1415f99a636bSHong Zhang } 1416f99a636bSHong Zhang 1417a3fe58edSHong Zhang ierr = MatDenseRestoreArray(Cden,&ca_den);CHKERRQ(ierr); 1418f99a636bSHong Zhang #if defined(PETSC_USE_INFO) 1419077f23c2SHong Zhang if (matcoloring->brows > 0) { 1420f99a636bSHong Zhang ierr = PetscInfo1(Csp,"Loop over %D row blocks for den2sp\n",brows);CHKERRQ(ierr); 1421e88777f2SHong Zhang } else { 1422077f23c2SHong Zhang ierr = PetscInfo(Csp,"Loop over colors/columns of Cden, inefficient for large sparse matrix product \n");CHKERRQ(ierr); 1423e88777f2SHong Zhang } 1424f99a636bSHong Zhang #endif 14258972f759SHong Zhang PetscFunctionReturn(0); 14268972f759SHong Zhang } 14278972f759SHong Zhang 1428b9af6bddSHong Zhang PetscErrorCode MatTransposeColoringCreate_SeqAIJ(Mat mat,ISColoring iscoloring,MatTransposeColoring c) 1429b1683b59SHong Zhang { 1430b1683b59SHong Zhang PetscErrorCode ierr; 1431e88777f2SHong Zhang PetscInt i,n,nrows,Nbs,j,k,m,ncols,col,cm; 14321a83f524SJed Brown const PetscInt *is,*ci,*cj,*row_idx; 1433b28d80bdSHong Zhang PetscInt nis = iscoloring->n,*rowhit,bs = 1; 1434b1683b59SHong Zhang IS *isa; 1435b28d80bdSHong Zhang Mat_SeqAIJ *csp = (Mat_SeqAIJ*)mat->data; 1436e88777f2SHong Zhang PetscInt *colorforrow,*rows,*rows_i,*idxhit,*spidx,*den2sp,*den2sp_i; 1437e88777f2SHong Zhang PetscInt *colorforcol,*columns,*columns_i,brows; 1438e88777f2SHong Zhang PetscBool flg; 1439b1683b59SHong Zhang 1440b1683b59SHong Zhang PetscFunctionBegin; 1441b1683b59SHong Zhang ierr = ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);CHKERRQ(ierr); 1442e99be685SHong Zhang 14437ecccc15SHong Zhang /* bs >1 is not being tested yet! */ 1444e88777f2SHong Zhang Nbs = mat->cmap->N/bs; 1445b1683b59SHong Zhang c->M = mat->rmap->N/bs; /* set total rows, columns and local rows */ 1446e88777f2SHong Zhang c->N = Nbs; 1447e88777f2SHong Zhang c->m = c->M; 1448b1683b59SHong Zhang c->rstart = 0; 1449e88777f2SHong Zhang c->brows = 100; 1450b1683b59SHong Zhang 1451b1683b59SHong Zhang c->ncolors = nis; 1452dcca6d9dSJed Brown ierr = PetscMalloc3(nis,&c->ncolumns,nis,&c->nrows,nis+1,&colorforrow);CHKERRQ(ierr); 1453854ce69bSBarry Smith ierr = PetscMalloc1(csp->nz+1,&rows);CHKERRQ(ierr); 1454854ce69bSBarry Smith ierr = PetscMalloc1(csp->nz+1,&den2sp);CHKERRQ(ierr); 1455e88777f2SHong Zhang 1456e88777f2SHong Zhang brows = c->brows; 1457c5929fdfSBarry Smith ierr = PetscOptionsGetInt(NULL,NULL,"-matden2sp_brows",&brows,&flg);CHKERRQ(ierr); 1458e88777f2SHong Zhang if (flg) c->brows = brows; 1459eeb4fd02SHong Zhang if (brows > 0) { 1460854ce69bSBarry Smith ierr = PetscMalloc1(nis+1,&c->lstart);CHKERRQ(ierr); 1461e88777f2SHong Zhang } 14622205254eSKarl Rupp 1463d2853540SHong Zhang colorforrow[0] = 0; 1464d2853540SHong Zhang rows_i = rows; 1465f99a636bSHong Zhang den2sp_i = den2sp; 1466b1683b59SHong Zhang 1467854ce69bSBarry Smith ierr = PetscMalloc1(nis+1,&colorforcol);CHKERRQ(ierr); 1468854ce69bSBarry Smith ierr = PetscMalloc1(Nbs+1,&columns);CHKERRQ(ierr); 14692205254eSKarl Rupp 1470d2853540SHong Zhang colorforcol[0] = 0; 1471d2853540SHong Zhang columns_i = columns; 1472d2853540SHong Zhang 1473077f23c2SHong Zhang /* get column-wise storage of mat */ 1474077f23c2SHong Zhang ierr = MatGetColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr); 1475b1683b59SHong Zhang 1476b28d80bdSHong Zhang cm = c->m; 1477854ce69bSBarry Smith ierr = PetscMalloc1(cm+1,&rowhit);CHKERRQ(ierr); 1478854ce69bSBarry Smith ierr = PetscMalloc1(cm+1,&idxhit);CHKERRQ(ierr); 1479f99a636bSHong Zhang for (i=0; i<nis; i++) { /* loop over color */ 1480b1683b59SHong Zhang ierr = ISGetLocalSize(isa[i],&n);CHKERRQ(ierr); 1481b1683b59SHong Zhang ierr = ISGetIndices(isa[i],&is);CHKERRQ(ierr); 14822205254eSKarl Rupp 1483b1683b59SHong Zhang c->ncolumns[i] = n; 1484b1683b59SHong Zhang if (n) { 1485d2853540SHong Zhang ierr = PetscMemcpy(columns_i,is,n*sizeof(PetscInt));CHKERRQ(ierr); 1486b1683b59SHong Zhang } 1487d2853540SHong Zhang colorforcol[i+1] = colorforcol[i] + n; 1488d2853540SHong Zhang columns_i += n; 1489b1683b59SHong Zhang 1490b1683b59SHong Zhang /* fast, crude version requires O(N*N) work */ 1491e8354b3eSHong Zhang ierr = PetscMemzero(rowhit,cm*sizeof(PetscInt));CHKERRQ(ierr); 1492e99be685SHong Zhang 1493b7caf3d6SHong Zhang for (j=0; j<n; j++) { /* loop over columns*/ 1494b1683b59SHong Zhang col = is[j]; 1495d2853540SHong Zhang row_idx = cj + ci[col]; 1496b1683b59SHong Zhang m = ci[col+1] - ci[col]; 1497b7caf3d6SHong Zhang for (k=0; k<m; k++) { /* loop over columns marking them in rowhit */ 1498e99be685SHong Zhang idxhit[*row_idx] = spidx[ci[col] + k]; 1499d2853540SHong Zhang rowhit[*row_idx++] = col + 1; 1500b1683b59SHong Zhang } 1501b1683b59SHong Zhang } 1502b1683b59SHong Zhang /* count the number of hits */ 1503b1683b59SHong Zhang nrows = 0; 1504e8354b3eSHong Zhang for (j=0; j<cm; j++) { 1505b1683b59SHong Zhang if (rowhit[j]) nrows++; 1506b1683b59SHong Zhang } 1507b1683b59SHong Zhang c->nrows[i] = nrows; 1508d2853540SHong Zhang colorforrow[i+1] = colorforrow[i] + nrows; 1509d2853540SHong Zhang 1510b1683b59SHong Zhang nrows = 0; 1511b7caf3d6SHong Zhang for (j=0; j<cm; j++) { /* loop over rows */ 1512b1683b59SHong Zhang if (rowhit[j]) { 1513d2853540SHong Zhang rows_i[nrows] = j; 151412b89a43SHong Zhang den2sp_i[nrows] = idxhit[j]; 1515b1683b59SHong Zhang nrows++; 1516b1683b59SHong Zhang } 1517b1683b59SHong Zhang } 1518e88777f2SHong Zhang den2sp_i += nrows; 1519077f23c2SHong Zhang 1520b1683b59SHong Zhang ierr = ISRestoreIndices(isa[i],&is);CHKERRQ(ierr); 1521f99a636bSHong Zhang rows_i += nrows; 1522b1683b59SHong Zhang } 15230298fd71SBarry Smith ierr = MatRestoreColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr); 1524b28d80bdSHong Zhang ierr = PetscFree(rowhit);CHKERRQ(ierr); 1525b1683b59SHong Zhang ierr = ISColoringRestoreIS(iscoloring,&isa);CHKERRQ(ierr); 1526d2853540SHong Zhang if (csp->nz != colorforrow[nis]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"csp->nz %d != colorforrow[nis] %d",csp->nz,colorforrow[nis]); 1527b28d80bdSHong Zhang 1528d2853540SHong Zhang c->colorforrow = colorforrow; 1529d2853540SHong Zhang c->rows = rows; 1530f99a636bSHong Zhang c->den2sp = den2sp; 1531d2853540SHong Zhang c->colorforcol = colorforcol; 1532d2853540SHong Zhang c->columns = columns; 15332205254eSKarl Rupp 1534f94ccd6cSHong Zhang ierr = PetscFree(idxhit);CHKERRQ(ierr); 1535b1683b59SHong Zhang PetscFunctionReturn(0); 1536b1683b59SHong Zhang } 1537d013fc79Sandi selinger 1538d013fc79Sandi selinger /* Needed for MatMatMult_SeqAIJ_SeqAIJ_Combined() */ 1539d013fc79Sandi selinger /* Append value to an array if the value is not present yet. A bitarray */ 1540d013fc79Sandi selinger /* was used to determine if there is already an entry at this position. */ 1541d013fc79Sandi selinger void appendToArray(PetscInt val, PetscInt *array, PetscInt *cnzi) 1542d013fc79Sandi selinger { 1543d013fc79Sandi selinger array[(*cnzi)++] = val; 1544d013fc79Sandi selinger } 1545d013fc79Sandi selinger 1546d013fc79Sandi selinger /* Needed for qsort in MatMatMult_SeqAIJ_SeqAIJ_Combined() */ 1547d013fc79Sandi selinger PetscInt cmpfunc (const void * a, const void * b) 1548d013fc79Sandi selinger { 1549d013fc79Sandi selinger return ( *(PetscInt*)a - *(PetscInt*)b ); 1550d013fc79Sandi selinger } 1551d013fc79Sandi selinger 1552*73b67375Sandi selinger /* This algorithm combines the symbolic and numeric phase of matrix-matrix multiplication. */ 1553d013fc79Sandi selinger PetscErrorCode MatMatMult_SeqAIJ_SeqAIJ_Combined(Mat A,Mat B,PetscReal fill,Mat *C) 1554d013fc79Sandi selinger { 1555d013fc79Sandi selinger PetscErrorCode ierr; 1556d013fc79Sandi selinger PetscLogDouble flops=0.0; 1557d013fc79Sandi selinger Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data, *b = (Mat_SeqAIJ*)B->data, *c; 1558d013fc79Sandi selinger const PetscInt *ai = a->i,*bi = b->i, *aj = a->j; 1559d013fc79Sandi selinger PetscInt *ci,*cj,*cj_i; 1560d013fc79Sandi selinger PetscScalar *ca, *ca_i; 1561d013fc79Sandi selinger PetscInt c_maxmem = 0, a_maxrownnz = 0, a_rownnz, a_col; 1562d013fc79Sandi selinger PetscInt am = A->rmap->N, bn = B->cmap->N, bm = B->rmap->N; 1563d013fc79Sandi selinger PetscInt i, k, ndouble = 0; 1564d013fc79Sandi selinger PetscReal afill; 1565d013fc79Sandi selinger PetscScalar *c_row_val_dense; 1566d013fc79Sandi selinger PetscBool *c_row_idx_flags; 1567d013fc79Sandi selinger PetscInt *aj_i = a->j; 1568d013fc79Sandi selinger PetscScalar *aa_i = a->a; 1569d013fc79Sandi selinger 1570d013fc79Sandi selinger PetscFunctionBegin; 1571d013fc79Sandi selinger /* Step 1: Determine upper bounds on memory for C */ 1572*73b67375Sandi selinger for (i=0; i<am; i++) { /* iterate over all rows of A */ 1573d013fc79Sandi selinger 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 */ 1574d013fc79Sandi selinger const PetscInt *acol = aj + ai[i]; /* column indices of nonzero entries in this row */ 1575d013fc79Sandi selinger a_rownnz = 0; 1576d013fc79Sandi selinger for (k=0;k<anzi;++k) a_rownnz += bi[acol[k]+1] - bi[acol[k]]; 1577d013fc79Sandi selinger a_maxrownnz = PetscMax(a_maxrownnz, a_rownnz); 1578d013fc79Sandi selinger c_maxmem += a_rownnz; 1579d013fc79Sandi selinger } 1580d013fc79Sandi selinger ierr = PetscMalloc1(am+1, &ci); CHKERRQ(ierr); 1581d013fc79Sandi selinger ierr = PetscMalloc1(bn, &c_row_val_dense); CHKERRQ(ierr); 1582d013fc79Sandi selinger ierr = PetscMalloc1(bn, &c_row_idx_flags); CHKERRQ(ierr); 1583d013fc79Sandi selinger ierr = PetscMalloc1(c_maxmem,&cj); CHKERRQ(ierr); 1584d013fc79Sandi selinger ierr = PetscMalloc1(c_maxmem,&ca); CHKERRQ(ierr); 1585d013fc79Sandi selinger ca_i = ca; 1586d013fc79Sandi selinger cj_i = cj; 1587d013fc79Sandi selinger ci[0] = 0; 1588*73b67375Sandi selinger ierr = PetscMemzero(c_row_val_dense, bn * sizeof(PetscScalar));CHKERRQ(ierr); 1589*73b67375Sandi selinger ierr = PetscMemzero(c_row_idx_flags, bn * sizeof(PetscBool));CHKERRQ(ierr); 1590d013fc79Sandi selinger for (i=0; i<am; i++) { 1591d013fc79Sandi selinger /* Step 2: Initialize the dense row vector for C */ 1592d013fc79Sandi selinger 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 */ 1593d013fc79Sandi selinger PetscInt cnzi = 0; 1594d013fc79Sandi selinger PetscInt *bj_i; 1595d013fc79Sandi selinger PetscScalar *ba_i; 1596d013fc79Sandi selinger 1597d013fc79Sandi selinger /* Step 3: Do the numerical calculations */ 1598d013fc79Sandi selinger for (a_col=0; a_col<anzi; a_col++) { /* iterate over all non zero values in a row of A */ 1599d013fc79Sandi selinger PetscInt a_col_index = aj_i[a_col]; 1600d013fc79Sandi selinger const PetscInt bnzi = bi[a_col_index+1] - bi[a_col_index]; 1601d013fc79Sandi selinger flops += 2*bnzi; 1602d013fc79Sandi selinger bj_i = b->j + bi[a_col_index]; /* points to the current row in bj */ 1603d013fc79Sandi selinger ba_i = b->a + bi[a_col_index]; /* points to the current row in ba */ 1604d013fc79Sandi selinger for (k=0; k<bnzi; ++k) { /* iterate over all non zeros of this row in B */ 1605d013fc79Sandi selinger if (c_row_idx_flags[ bj_i[k] ] == PETSC_FALSE) { 1606d013fc79Sandi selinger appendToArray(bj_i[k], cj_i, &cnzi); 1607d013fc79Sandi selinger c_row_idx_flags[ bj_i[k] ] = PETSC_TRUE; 1608d013fc79Sandi selinger } 1609d013fc79Sandi selinger c_row_val_dense[ bj_i[k] ] += aa_i[a_col] * ba_i[k]; 1610d013fc79Sandi selinger } 1611d013fc79Sandi selinger } 1612d013fc79Sandi selinger 1613d013fc79Sandi selinger /* Sort array */ 1614d013fc79Sandi selinger qsort(cj_i, cnzi, sizeof(PetscInt), cmpfunc); 1615d013fc79Sandi selinger /* Step 4 */ 1616d013fc79Sandi selinger for (k=0; k < cnzi; k++) { 1617d013fc79Sandi selinger ca_i[k] = c_row_val_dense[cj_i[k]]; 1618d013fc79Sandi selinger c_row_val_dense[cj_i[k]] = 0.; 1619d013fc79Sandi selinger c_row_idx_flags[cj_i[k]] = PETSC_FALSE; 1620d013fc79Sandi selinger } 1621d013fc79Sandi selinger /* terminate current row */ 1622d013fc79Sandi selinger aa_i += anzi; 1623d013fc79Sandi selinger aj_i += anzi; 1624d013fc79Sandi selinger ca_i += cnzi; 1625d013fc79Sandi selinger cj_i += cnzi; 1626d013fc79Sandi selinger ci[i+1] = ci[i] + cnzi; 1627d013fc79Sandi selinger flops += cnzi; 1628d013fc79Sandi selinger } 1629d013fc79Sandi selinger 1630d013fc79Sandi selinger /* Step 5 */ 1631d013fc79Sandi selinger /* Create the new matrix */ 1632d013fc79Sandi selinger ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,NULL,C);CHKERRQ(ierr); 1633d013fc79Sandi selinger ierr = MatSetBlockSizesFromMats(*C,A,B);CHKERRQ(ierr); 1634d013fc79Sandi selinger 1635d013fc79Sandi selinger /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 1636d013fc79Sandi selinger /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */ 1637d013fc79Sandi selinger c = (Mat_SeqAIJ*)((*C)->data); 1638d013fc79Sandi selinger c->a = ca; 1639d013fc79Sandi selinger c->free_a = PETSC_TRUE; 1640d013fc79Sandi selinger c->free_ij = PETSC_TRUE; 1641d013fc79Sandi selinger c->nonew = 0; 1642d013fc79Sandi selinger 1643d013fc79Sandi selinger /* set MatInfo */ 1644d013fc79Sandi selinger afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5; 1645d013fc79Sandi selinger if (afill < 1.0) afill = 1.0; 1646d013fc79Sandi selinger c->maxnz = ci[am]; 1647d013fc79Sandi selinger c->nz = ci[am]; 1648d013fc79Sandi selinger (*C)->info.mallocs = ndouble; 1649d013fc79Sandi selinger (*C)->info.fill_ratio_given = fill; 1650d013fc79Sandi selinger (*C)->info.fill_ratio_needed = afill; 1651d013fc79Sandi selinger 1652*73b67375Sandi selinger ierr = PetscFree(c_row_val_dense);CHKERRQ(ierr); 1653*73b67375Sandi selinger ierr = PetscFree(c_row_idx_flags);CHKERRQ(ierr); 1654d013fc79Sandi selinger 1655d013fc79Sandi selinger ierr = MatAssemblyBegin(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1656d013fc79Sandi selinger ierr = MatAssemblyEnd(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1657d013fc79Sandi selinger ierr = PetscLogFlops(flops);CHKERRQ(ierr); 1658d013fc79Sandi selinger PetscFunctionReturn(0); 1659d013fc79Sandi selinger } 1660