xref: /petsc/src/mat/impls/aij/seq/matmatmult.c (revision 02fe19650ff8cf45c2e973b4ee0f4cdf7bff9c74)
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),&current_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);
152*02fe1965SBarry Smith    ierr = MatSetType(*C,((PetscObject)A)->type_name);CHKERRQ(ierr);
15358c24d83SHong Zhang 
15458c24d83SHong Zhang   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
15558c24d83SHong Zhang   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
15658c24d83SHong Zhang   c                         = (Mat_SeqAIJ*)((*C)->data);
157aa1aec99SHong Zhang   c->free_a                 = PETSC_FALSE;
158e6b907acSBarry Smith   c->free_ij                = PETSC_TRUE;
15958c24d83SHong Zhang   c->nonew                  = 0;
160002e8affSHong Zhang   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ; /* fast, needs non-scalable O(bn) array 'abdense' */
1610b7e3e3dSHong Zhang 
1627212da7cSHong Zhang   /* set MatInfo */
1637212da7cSHong Zhang   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
164f2b054eeSHong Zhang   if (afill < 1.0) afill = 1.0;
1657212da7cSHong Zhang   c->maxnz                     = ci[am];
1667212da7cSHong Zhang   c->nz                        = ci[am];
167fb908031SHong Zhang   (*C)->info.mallocs           = ndouble;
1687212da7cSHong Zhang   (*C)->info.fill_ratio_given  = fill;
1697212da7cSHong Zhang   (*C)->info.fill_ratio_needed = afill;
1707212da7cSHong Zhang 
1717212da7cSHong Zhang #if defined(PETSC_USE_INFO)
1727212da7cSHong Zhang   if (ci[am]) {
17357622a8eSBarry Smith     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %g needed %g.\n",ndouble,(double)fill,(double)afill);CHKERRQ(ierr);
17457622a8eSBarry Smith     ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr);
175f2b054eeSHong Zhang   } else {
176f2b054eeSHong Zhang     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
177be0fcf8dSHong Zhang   }
178f2b054eeSHong Zhang #endif
17958c24d83SHong Zhang   PetscFunctionReturn(0);
18058c24d83SHong Zhang }
181d50806bdSBarry Smith 
182dfbe8321SBarry Smith PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
183d50806bdSBarry Smith {
184dfbe8321SBarry Smith   PetscErrorCode ierr;
185d13dce4bSSatish Balay   PetscLogDouble flops=0.0;
186d50806bdSBarry Smith   Mat_SeqAIJ     *a   = (Mat_SeqAIJ*)A->data;
187d50806bdSBarry Smith   Mat_SeqAIJ     *b   = (Mat_SeqAIJ*)B->data;
188d50806bdSBarry Smith   Mat_SeqAIJ     *c   = (Mat_SeqAIJ*)C->data;
18938baddfdSBarry Smith   PetscInt       *ai  =a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j;
190c58ca2e3SHong Zhang   PetscInt       am   =A->rmap->n,cm=C->rmap->n;
191519eb980SHong Zhang   PetscInt       i,j,k,anzi,bnzi,cnzi,brow;
192aa1aec99SHong Zhang   PetscScalar    *aa=a->a,*ba=b->a,*baj,*ca,valtmp;
193aa1aec99SHong Zhang   PetscScalar    *ab_dense;
194d50806bdSBarry Smith 
195d50806bdSBarry Smith   PetscFunctionBegin;
196aa1aec99SHong Zhang   if (!c->a) { /* first call of MatMatMultNumeric_SeqAIJ_SeqAIJ, allocate ca and matmult_abdense */
197854ce69bSBarry Smith     ierr      = PetscMalloc1(ci[cm]+1,&ca);CHKERRQ(ierr);
198aa1aec99SHong Zhang     c->a      = ca;
199aa1aec99SHong Zhang     c->free_a = PETSC_TRUE;
200aa1aec99SHong Zhang   } else {
201aa1aec99SHong Zhang     ca        = c->a;
202d90d86d0SMatthew G. Knepley   }
203d90d86d0SMatthew G. Knepley   if (!c->matmult_abdense) {
2041795a4d1SJed Brown     ierr = PetscCalloc1(B->cmap->N,&ab_dense);CHKERRQ(ierr);
205d90d86d0SMatthew G. Knepley     c->matmult_abdense = ab_dense;
206d90d86d0SMatthew G. Knepley   } else {
207aa1aec99SHong Zhang     ab_dense = c->matmult_abdense;
208aa1aec99SHong Zhang   }
209aa1aec99SHong Zhang 
210c124e916SHong Zhang   /* clean old values in C */
211c124e916SHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
212d50806bdSBarry Smith   /* Traverse A row-wise. */
213d50806bdSBarry Smith   /* Build the ith row in C by summing over nonzero columns in A, */
214d50806bdSBarry Smith   /* the rows of B corresponding to nonzeros of A. */
215d50806bdSBarry Smith   for (i=0; i<am; i++) {
216d50806bdSBarry Smith     anzi = ai[i+1] - ai[i];
217d50806bdSBarry Smith     for (j=0; j<anzi; j++) {
218519eb980SHong Zhang       brow = aj[j];
219d50806bdSBarry Smith       bnzi = bi[brow+1] - bi[brow];
220d50806bdSBarry Smith       bjj  = bj + bi[brow];
221d50806bdSBarry Smith       baj  = ba + bi[brow];
222519eb980SHong Zhang       /* perform dense axpy */
22336ec6d2dSHong Zhang       valtmp = aa[j];
224519eb980SHong Zhang       for (k=0; k<bnzi; k++) {
22536ec6d2dSHong Zhang         ab_dense[bjj[k]] += valtmp*baj[k];
226519eb980SHong Zhang       }
227519eb980SHong Zhang       flops += 2*bnzi;
228519eb980SHong Zhang     }
229c58ca2e3SHong Zhang     aj += anzi; aa += anzi;
230c58ca2e3SHong Zhang 
231c58ca2e3SHong Zhang     cnzi = ci[i+1] - ci[i];
232519eb980SHong Zhang     for (k=0; k<cnzi; k++) {
233519eb980SHong Zhang       ca[k]          += ab_dense[cj[k]];
234519eb980SHong Zhang       ab_dense[cj[k]] = 0.0; /* zero ab_dense */
235519eb980SHong Zhang     }
236519eb980SHong Zhang     flops += cnzi;
237519eb980SHong Zhang     cj    += cnzi; ca += cnzi;
238519eb980SHong Zhang   }
239c58ca2e3SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
240c58ca2e3SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
241c58ca2e3SHong Zhang   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
242c58ca2e3SHong Zhang   PetscFunctionReturn(0);
243c58ca2e3SHong Zhang }
244c58ca2e3SHong Zhang 
24525023028SHong Zhang PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(Mat A,Mat B,Mat C)
246c58ca2e3SHong Zhang {
247c58ca2e3SHong Zhang   PetscErrorCode ierr;
248c58ca2e3SHong Zhang   PetscLogDouble flops=0.0;
249c58ca2e3SHong Zhang   Mat_SeqAIJ     *a   = (Mat_SeqAIJ*)A->data;
250c58ca2e3SHong Zhang   Mat_SeqAIJ     *b   = (Mat_SeqAIJ*)B->data;
251c58ca2e3SHong Zhang   Mat_SeqAIJ     *c   = (Mat_SeqAIJ*)C->data;
252c58ca2e3SHong Zhang   PetscInt       *ai  = a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j;
253c58ca2e3SHong Zhang   PetscInt       am   = A->rmap->N,cm=C->rmap->N;
254c58ca2e3SHong Zhang   PetscInt       i,j,k,anzi,bnzi,cnzi,brow;
25536ec6d2dSHong Zhang   PetscScalar    *aa=a->a,*ba=b->a,*baj,*ca=c->a,valtmp;
256c58ca2e3SHong Zhang   PetscInt       nextb;
257c58ca2e3SHong Zhang 
258c58ca2e3SHong Zhang   PetscFunctionBegin;
259cf742fccSHong Zhang   if (!ca) { /* first call of MatMatMultNumeric_SeqAIJ_SeqAIJ, allocate ca and matmult_abdense */
260cf742fccSHong Zhang     ierr      = PetscMalloc1(ci[cm]+1,&ca);CHKERRQ(ierr);
261cf742fccSHong Zhang     c->a      = ca;
262cf742fccSHong Zhang     c->free_a = PETSC_TRUE;
263cf742fccSHong Zhang   }
264cf742fccSHong Zhang 
265c58ca2e3SHong Zhang   /* clean old values in C */
266c58ca2e3SHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
267c58ca2e3SHong Zhang   /* Traverse A row-wise. */
268c58ca2e3SHong Zhang   /* Build the ith row in C by summing over nonzero columns in A, */
269c58ca2e3SHong Zhang   /* the rows of B corresponding to nonzeros of A. */
270519eb980SHong Zhang   for (i=0; i<am; i++) {
271519eb980SHong Zhang     anzi = ai[i+1] - ai[i];
272519eb980SHong Zhang     cnzi = ci[i+1] - ci[i];
273519eb980SHong Zhang     for (j=0; j<anzi; j++) {
274519eb980SHong Zhang       brow = aj[j];
275519eb980SHong Zhang       bnzi = bi[brow+1] - bi[brow];
276519eb980SHong Zhang       bjj  = bj + bi[brow];
277519eb980SHong Zhang       baj  = ba + bi[brow];
278519eb980SHong Zhang       /* perform sparse axpy */
27936ec6d2dSHong Zhang       valtmp = aa[j];
280c124e916SHong Zhang       nextb  = 0;
281c124e916SHong Zhang       for (k=0; nextb<bnzi; k++) {
282c124e916SHong Zhang         if (cj[k] == bjj[nextb]) { /* ccol == bcol */
28336ec6d2dSHong Zhang           ca[k] += valtmp*baj[nextb++];
284c124e916SHong Zhang         }
285d50806bdSBarry Smith       }
286d50806bdSBarry Smith       flops += 2*bnzi;
287d50806bdSBarry Smith     }
288519eb980SHong Zhang     aj += anzi; aa += anzi;
289519eb980SHong Zhang     cj += cnzi; ca += cnzi;
290d50806bdSBarry Smith   }
291c58ca2e3SHong Zhang 
292716bacf3SKris Buschelman   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
293716bacf3SKris Buschelman   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
294d50806bdSBarry Smith   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
295d50806bdSBarry Smith   PetscFunctionReturn(0);
296d50806bdSBarry Smith }
297bc011b1eSHong Zhang 
2983c50cad2SHong Zhang PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast(Mat A,Mat B,PetscReal fill,Mat *C)
29925296bd5SBarry Smith {
30025296bd5SBarry Smith   PetscErrorCode     ierr;
30125296bd5SBarry Smith   Mat_SeqAIJ         *a  = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
30225296bd5SBarry Smith   PetscInt           *ai = a->i,*bi=b->i,*ci,*cj;
3033c50cad2SHong Zhang   PetscInt           am  = A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
30425296bd5SBarry Smith   MatScalar          *ca;
30525296bd5SBarry Smith   PetscReal          afill;
306eca6b86aSHong Zhang   PetscInt           i,j,anzi,brow,bnzj,cnzi,*bj,*aj,*lnk,ndouble=0,Crmax;
307eca6b86aSHong Zhang   PetscTable         ta;
3080298fd71SBarry Smith   PetscFreeSpaceList free_space=NULL,current_space=NULL;
30925296bd5SBarry Smith 
31025296bd5SBarry Smith   PetscFunctionBegin;
3113c50cad2SHong Zhang   /* Get ci and cj - same as MatMatMultSymbolic_SeqAIJ_SeqAIJ except using PetscLLxxx_fast() */
3123c50cad2SHong Zhang   /*-----------------------------------------------------------------------------------------*/
3133c50cad2SHong Zhang   /* Allocate arrays for fill computation and free space for accumulating nonzero column */
314854ce69bSBarry Smith   ierr  = PetscMalloc1(am+2,&ci);CHKERRQ(ierr);
3153c50cad2SHong Zhang   ci[0] = 0;
3163c50cad2SHong Zhang 
3173c50cad2SHong Zhang   /* create and initialize a linked list */
318c373ccc6SHong Zhang   ierr = PetscTableCreate(bn,bn,&ta);CHKERRQ(ierr);
319c373ccc6SHong Zhang   MatRowMergeMax_SeqAIJ(b,bm,ta);
320eca6b86aSHong Zhang   ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr);
321eca6b86aSHong Zhang   ierr = PetscTableDestroy(&ta);CHKERRQ(ierr);
322eca6b86aSHong Zhang 
323eca6b86aSHong Zhang   ierr = PetscLLCondensedCreate_fast(Crmax,&lnk);CHKERRQ(ierr);
3243c50cad2SHong Zhang 
3253c50cad2SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
326f91af8c7SBarry Smith   ierr          = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ai[am],bi[bm])),&free_space);CHKERRQ(ierr);
3273c50cad2SHong Zhang   current_space = free_space;
3283c50cad2SHong Zhang 
3293c50cad2SHong Zhang   /* Determine ci and cj */
3303c50cad2SHong Zhang   for (i=0; i<am; i++) {
3313c50cad2SHong Zhang     anzi = ai[i+1] - ai[i];
3323c50cad2SHong Zhang     aj   = a->j + ai[i];
3333c50cad2SHong Zhang     for (j=0; j<anzi; j++) {
3343c50cad2SHong Zhang       brow = aj[j];
3353c50cad2SHong Zhang       bnzj = bi[brow+1] - bi[brow];
3363c50cad2SHong Zhang       bj   = b->j + bi[brow];
3373c50cad2SHong Zhang       /* add non-zero cols of B into the sorted linked list lnk */
3383c50cad2SHong Zhang       ierr = PetscLLCondensedAddSorted_fast(bnzj,bj,lnk);CHKERRQ(ierr);
3393c50cad2SHong Zhang     }
3403c50cad2SHong Zhang     cnzi = lnk[1];
3413c50cad2SHong Zhang 
3423c50cad2SHong Zhang     /* If free space is not available, make more free space */
3433c50cad2SHong Zhang     /* Double the amount of total space in the list */
3443c50cad2SHong Zhang     if (current_space->local_remaining<cnzi) {
345f91af8c7SBarry Smith       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(cnzi,current_space->total_array_size),&current_space);CHKERRQ(ierr);
3463c50cad2SHong Zhang       ndouble++;
3473c50cad2SHong Zhang     }
3483c50cad2SHong Zhang 
3493c50cad2SHong Zhang     /* Copy data into free space, then initialize lnk */
3503c50cad2SHong Zhang     ierr = PetscLLCondensedClean_fast(cnzi,current_space->array,lnk);CHKERRQ(ierr);
3512205254eSKarl Rupp 
3523c50cad2SHong Zhang     current_space->array           += cnzi;
3533c50cad2SHong Zhang     current_space->local_used      += cnzi;
3543c50cad2SHong Zhang     current_space->local_remaining -= cnzi;
3552205254eSKarl Rupp 
3563c50cad2SHong Zhang     ci[i+1] = ci[i] + cnzi;
3573c50cad2SHong Zhang   }
3583c50cad2SHong Zhang 
3593c50cad2SHong Zhang   /* Column indices are in the list of free space */
3603c50cad2SHong Zhang   /* Allocate space for cj, initialize cj, and */
3613c50cad2SHong Zhang   /* destroy list of free space and other temporary array(s) */
362854ce69bSBarry Smith   ierr = PetscMalloc1(ci[am]+1,&cj);CHKERRQ(ierr);
3633c50cad2SHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
3643c50cad2SHong Zhang   ierr = PetscLLCondensedDestroy_fast(lnk);CHKERRQ(ierr);
36525296bd5SBarry Smith 
36625296bd5SBarry Smith   /* Allocate space for ca */
367854ce69bSBarry Smith   ierr = PetscMalloc1(ci[am]+1,&ca);CHKERRQ(ierr);
36825296bd5SBarry Smith   ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr);
36925296bd5SBarry Smith 
37025296bd5SBarry Smith   /* put together the new symbolic matrix */
371ce94432eSBarry Smith   ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,ca,C);CHKERRQ(ierr);
37233d57670SJed Brown   ierr = MatSetBlockSizesFromMats(*C,A,B);CHKERRQ(ierr);
373*02fe1965SBarry Smith   ierr = MatSetType(*C,((PetscObject)A)->type_name);CHKERRQ(ierr);
37425296bd5SBarry Smith 
37525296bd5SBarry Smith   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
37625296bd5SBarry Smith   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
37725296bd5SBarry Smith   c          = (Mat_SeqAIJ*)((*C)->data);
37825296bd5SBarry Smith   c->free_a  = PETSC_TRUE;
37925296bd5SBarry Smith   c->free_ij = PETSC_TRUE;
38025296bd5SBarry Smith   c->nonew   = 0;
3812205254eSKarl Rupp 
38225296bd5SBarry Smith   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable; /* slower, less memory */
38325296bd5SBarry Smith 
38425296bd5SBarry Smith   /* set MatInfo */
38525296bd5SBarry Smith   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
38625296bd5SBarry Smith   if (afill < 1.0) afill = 1.0;
38725296bd5SBarry Smith   c->maxnz                     = ci[am];
38825296bd5SBarry Smith   c->nz                        = ci[am];
3893c50cad2SHong Zhang   (*C)->info.mallocs           = ndouble;
39025296bd5SBarry Smith   (*C)->info.fill_ratio_given  = fill;
39125296bd5SBarry Smith   (*C)->info.fill_ratio_needed = afill;
39225296bd5SBarry Smith 
39325296bd5SBarry Smith #if defined(PETSC_USE_INFO)
39425296bd5SBarry Smith   if (ci[am]) {
39557622a8eSBarry Smith     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %g needed %g.\n",ndouble,(double)fill,(double)afill);CHKERRQ(ierr);
39657622a8eSBarry Smith     ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr);
39725296bd5SBarry Smith   } else {
39825296bd5SBarry Smith     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
39925296bd5SBarry Smith   }
40025296bd5SBarry Smith #endif
40125296bd5SBarry Smith   PetscFunctionReturn(0);
40225296bd5SBarry Smith }
40325296bd5SBarry Smith 
40425296bd5SBarry Smith 
40525023028SHong Zhang PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(Mat A,Mat B,PetscReal fill,Mat *C)
406e9e4536cSHong Zhang {
407e9e4536cSHong Zhang   PetscErrorCode     ierr;
408e9e4536cSHong Zhang   Mat_SeqAIJ         *a  = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
409bf9555e6SHong Zhang   PetscInt           *ai = a->i,*bi=b->i,*ci,*cj;
41025c41797SHong Zhang   PetscInt           am  = A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
411e9e4536cSHong Zhang   MatScalar          *ca;
412e9e4536cSHong Zhang   PetscReal          afill;
413eca6b86aSHong Zhang   PetscInt           i,j,anzi,brow,bnzj,cnzi,*bj,*aj,*lnk,ndouble=0,Crmax;
414eca6b86aSHong Zhang   PetscTable         ta;
4150298fd71SBarry Smith   PetscFreeSpaceList free_space=NULL,current_space=NULL;
416e9e4536cSHong Zhang 
417e9e4536cSHong Zhang   PetscFunctionBegin;
418bd958071SHong Zhang   /* Get ci and cj - same as MatMatMultSymbolic_SeqAIJ_SeqAIJ except using PetscLLxxx_Scalalbe() */
419bd958071SHong Zhang   /*---------------------------------------------------------------------------------------------*/
420bd958071SHong Zhang   /* Allocate arrays for fill computation and free space for accumulating nonzero column */
421854ce69bSBarry Smith   ierr  = PetscMalloc1(am+2,&ci);CHKERRQ(ierr);
422bd958071SHong Zhang   ci[0] = 0;
423bd958071SHong Zhang 
424bd958071SHong Zhang   /* create and initialize a linked list */
425c373ccc6SHong Zhang   ierr = PetscTableCreate(bn,bn,&ta);CHKERRQ(ierr);
426c373ccc6SHong Zhang   MatRowMergeMax_SeqAIJ(b,bm,ta);
427eca6b86aSHong Zhang   ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr);
428eca6b86aSHong Zhang   ierr = PetscTableDestroy(&ta);CHKERRQ(ierr);
429eca6b86aSHong Zhang   ierr = PetscLLCondensedCreate_Scalable(Crmax,&lnk);CHKERRQ(ierr);
430bd958071SHong Zhang 
431bd958071SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
432f91af8c7SBarry Smith   ierr          = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ai[am],bi[bm])),&free_space);CHKERRQ(ierr);
433bd958071SHong Zhang   current_space = free_space;
434bd958071SHong Zhang 
435bd958071SHong Zhang   /* Determine ci and cj */
436bd958071SHong Zhang   for (i=0; i<am; i++) {
437bd958071SHong Zhang     anzi = ai[i+1] - ai[i];
438bd958071SHong Zhang     aj   = a->j + ai[i];
439bd958071SHong Zhang     for (j=0; j<anzi; j++) {
440bd958071SHong Zhang       brow = aj[j];
441bd958071SHong Zhang       bnzj = bi[brow+1] - bi[brow];
442bd958071SHong Zhang       bj   = b->j + bi[brow];
443bd958071SHong Zhang       /* add non-zero cols of B into the sorted linked list lnk */
444bd958071SHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(bnzj,bj,lnk);CHKERRQ(ierr);
445bd958071SHong Zhang     }
446bd958071SHong Zhang     cnzi = lnk[0];
447bd958071SHong Zhang 
448bd958071SHong Zhang     /* If free space is not available, make more free space */
449bd958071SHong Zhang     /* Double the amount of total space in the list */
450bd958071SHong Zhang     if (current_space->local_remaining<cnzi) {
451f91af8c7SBarry Smith       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(cnzi,current_space->total_array_size),&current_space);CHKERRQ(ierr);
452bd958071SHong Zhang       ndouble++;
453bd958071SHong Zhang     }
454bd958071SHong Zhang 
455bd958071SHong Zhang     /* Copy data into free space, then initialize lnk */
456bd958071SHong Zhang     ierr = PetscLLCondensedClean_Scalable(cnzi,current_space->array,lnk);CHKERRQ(ierr);
4572205254eSKarl Rupp 
458bd958071SHong Zhang     current_space->array           += cnzi;
459bd958071SHong Zhang     current_space->local_used      += cnzi;
460bd958071SHong Zhang     current_space->local_remaining -= cnzi;
4612205254eSKarl Rupp 
462bd958071SHong Zhang     ci[i+1] = ci[i] + cnzi;
463bd958071SHong Zhang   }
464bd958071SHong Zhang 
465bd958071SHong Zhang   /* Column indices are in the list of free space */
466bd958071SHong Zhang   /* Allocate space for cj, initialize cj, and */
467bd958071SHong Zhang   /* destroy list of free space and other temporary array(s) */
468854ce69bSBarry Smith   ierr = PetscMalloc1(ci[am]+1,&cj);CHKERRQ(ierr);
469bd958071SHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
470bd958071SHong Zhang   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr);
471e9e4536cSHong Zhang 
472e9e4536cSHong Zhang   /* Allocate space for ca */
473bd958071SHong Zhang   /*-----------------------*/
474854ce69bSBarry Smith   ierr = PetscMalloc1(ci[am]+1,&ca);CHKERRQ(ierr);
475e9e4536cSHong Zhang   ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr);
476e9e4536cSHong Zhang 
477e9e4536cSHong Zhang   /* put together the new symbolic matrix */
478ce94432eSBarry Smith   ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,ca,C);CHKERRQ(ierr);
47933d57670SJed Brown   ierr = MatSetBlockSizesFromMats(*C,A,B);CHKERRQ(ierr);
480*02fe1965SBarry Smith   ierr = MatSetType(*C,((PetscObject)A)->type_name);CHKERRQ(ierr);
481e9e4536cSHong Zhang 
482e9e4536cSHong Zhang   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
483e9e4536cSHong Zhang   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
484e9e4536cSHong Zhang   c          = (Mat_SeqAIJ*)((*C)->data);
485e9e4536cSHong Zhang   c->free_a  = PETSC_TRUE;
486e9e4536cSHong Zhang   c->free_ij = PETSC_TRUE;
487e9e4536cSHong Zhang   c->nonew   = 0;
4882205254eSKarl Rupp 
48925023028SHong Zhang   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable; /* slower, less memory */
490e9e4536cSHong Zhang 
491e9e4536cSHong Zhang   /* set MatInfo */
492e9e4536cSHong Zhang   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
493e9e4536cSHong Zhang   if (afill < 1.0) afill = 1.0;
494e9e4536cSHong Zhang   c->maxnz                     = ci[am];
495e9e4536cSHong Zhang   c->nz                        = ci[am];
496bd958071SHong Zhang   (*C)->info.mallocs           = ndouble;
497e9e4536cSHong Zhang   (*C)->info.fill_ratio_given  = fill;
498e9e4536cSHong Zhang   (*C)->info.fill_ratio_needed = afill;
499e9e4536cSHong Zhang 
500e9e4536cSHong Zhang #if defined(PETSC_USE_INFO)
501e9e4536cSHong Zhang   if (ci[am]) {
50257622a8eSBarry Smith     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %g needed %g.\n",ndouble,(double)fill,(double)afill);CHKERRQ(ierr);
50357622a8eSBarry Smith     ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr);
504e9e4536cSHong Zhang   } else {
505e9e4536cSHong Zhang     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
506e9e4536cSHong Zhang   }
507e9e4536cSHong Zhang #endif
508e9e4536cSHong Zhang   PetscFunctionReturn(0);
509e9e4536cSHong Zhang }
510e9e4536cSHong Zhang 
5110ced3a2bSJed Brown PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Heap(Mat A,Mat B,PetscReal fill,Mat *C)
5120ced3a2bSJed Brown {
5130ced3a2bSJed Brown   PetscErrorCode     ierr;
5140ced3a2bSJed Brown   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
5150ced3a2bSJed Brown   const PetscInt     *ai=a->i,*bi=b->i,*aj=a->j,*bj=b->j;
5160ced3a2bSJed Brown   PetscInt           *ci,*cj,*bb;
5170ced3a2bSJed Brown   PetscInt           am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
5180ced3a2bSJed Brown   PetscReal          afill;
5190ced3a2bSJed Brown   PetscInt           i,j,col,ndouble = 0;
5200298fd71SBarry Smith   PetscFreeSpaceList free_space=NULL,current_space=NULL;
5210ced3a2bSJed Brown   PetscHeap          h;
5220ced3a2bSJed Brown 
5230ced3a2bSJed Brown   PetscFunctionBegin;
524cd093f1dSJed Brown   /* Get ci and cj - by merging sorted rows using a heap */
5250ced3a2bSJed Brown   /*---------------------------------------------------------------------------------------------*/
5260ced3a2bSJed Brown   /* Allocate arrays for fill computation and free space for accumulating nonzero column */
527854ce69bSBarry Smith   ierr  = PetscMalloc1(am+2,&ci);CHKERRQ(ierr);
5280ced3a2bSJed Brown   ci[0] = 0;
5290ced3a2bSJed Brown 
5300ced3a2bSJed Brown   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
531f91af8c7SBarry Smith   ierr          = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ai[am],bi[bm])),&free_space);CHKERRQ(ierr);
5320ced3a2bSJed Brown   current_space = free_space;
5330ced3a2bSJed Brown 
5340ced3a2bSJed Brown   ierr = PetscHeapCreate(a->rmax,&h);CHKERRQ(ierr);
535785e854fSJed Brown   ierr = PetscMalloc1(a->rmax,&bb);CHKERRQ(ierr);
5360ced3a2bSJed Brown 
5370ced3a2bSJed Brown   /* Determine ci and cj */
5380ced3a2bSJed Brown   for (i=0; i<am; i++) {
5390ced3a2bSJed 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 */
5400ced3a2bSJed Brown     const PetscInt *acol = aj + ai[i]; /* column indices of nonzero entries in this row */
5410ced3a2bSJed Brown     ci[i+1] = ci[i];
5420ced3a2bSJed Brown     /* Populate the min heap */
5430ced3a2bSJed Brown     for (j=0; j<anzi; j++) {
5440ced3a2bSJed Brown       bb[j] = bi[acol[j]];         /* bb points at the start of the row */
5450ced3a2bSJed Brown       if (bb[j] < bi[acol[j]+1]) { /* Add if row is nonempty */
5460ced3a2bSJed Brown         ierr = PetscHeapAdd(h,j,bj[bb[j]++]);CHKERRQ(ierr);
5470ced3a2bSJed Brown       }
5480ced3a2bSJed Brown     }
5490ced3a2bSJed Brown     /* Pick off the min element, adding it to free space */
5500ced3a2bSJed Brown     ierr = PetscHeapPop(h,&j,&col);CHKERRQ(ierr);
5510ced3a2bSJed Brown     while (j >= 0) {
5520ced3a2bSJed Brown       if (current_space->local_remaining < 1) { /* double the size, but don't exceed 16 MiB */
553f91af8c7SBarry Smith         ierr = PetscFreeSpaceGet(PetscMin(PetscIntMultTruncate(2,current_space->total_array_size),16 << 20),&current_space);CHKERRQ(ierr);
5540ced3a2bSJed Brown         ndouble++;
5550ced3a2bSJed Brown       }
5560ced3a2bSJed Brown       *(current_space->array++) = col;
5570ced3a2bSJed Brown       current_space->local_used++;
5580ced3a2bSJed Brown       current_space->local_remaining--;
5590ced3a2bSJed Brown       ci[i+1]++;
5600ced3a2bSJed Brown 
5610ced3a2bSJed Brown       /* stash if anything else remains in this row of B */
5620ced3a2bSJed Brown       if (bb[j] < bi[acol[j]+1]) {ierr = PetscHeapStash(h,j,bj[bb[j]++]);CHKERRQ(ierr);}
5630ced3a2bSJed Brown       while (1) {               /* pop and stash any other rows of B that also had an entry in this column */
5640ced3a2bSJed Brown         PetscInt j2,col2;
5650ced3a2bSJed Brown         ierr = PetscHeapPeek(h,&j2,&col2);CHKERRQ(ierr);
5660ced3a2bSJed Brown         if (col2 != col) break;
5670ced3a2bSJed Brown         ierr = PetscHeapPop(h,&j2,&col2);CHKERRQ(ierr);
5680ced3a2bSJed Brown         if (bb[j2] < bi[acol[j2]+1]) {ierr = PetscHeapStash(h,j2,bj[bb[j2]++]);CHKERRQ(ierr);}
5690ced3a2bSJed Brown       }
5700ced3a2bSJed Brown       /* Put any stashed elements back into the min heap */
5710ced3a2bSJed Brown       ierr = PetscHeapUnstash(h);CHKERRQ(ierr);
5720ced3a2bSJed Brown       ierr = PetscHeapPop(h,&j,&col);CHKERRQ(ierr);
5730ced3a2bSJed Brown     }
5740ced3a2bSJed Brown   }
5750ced3a2bSJed Brown   ierr = PetscFree(bb);CHKERRQ(ierr);
5760ced3a2bSJed Brown   ierr = PetscHeapDestroy(&h);CHKERRQ(ierr);
5770ced3a2bSJed Brown 
5780ced3a2bSJed Brown   /* Column indices are in the list of free space */
5790ced3a2bSJed Brown   /* Allocate space for cj, initialize cj, and */
5800ced3a2bSJed Brown   /* destroy list of free space and other temporary array(s) */
581785e854fSJed Brown   ierr = PetscMalloc1(ci[am],&cj);CHKERRQ(ierr);
5820ced3a2bSJed Brown   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
5830ced3a2bSJed Brown 
5840ced3a2bSJed Brown   /* put together the new symbolic matrix */
585ce94432eSBarry Smith   ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,NULL,C);CHKERRQ(ierr);
58633d57670SJed Brown   ierr = MatSetBlockSizesFromMats(*C,A,B);CHKERRQ(ierr);
587*02fe1965SBarry Smith   ierr = MatSetType(*C,((PetscObject)A)->type_name);CHKERRQ(ierr);
5880ced3a2bSJed Brown 
5890ced3a2bSJed Brown   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
5900ced3a2bSJed Brown   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
5910ced3a2bSJed Brown   c          = (Mat_SeqAIJ*)((*C)->data);
5920ced3a2bSJed Brown   c->free_a  = PETSC_TRUE;
5930ced3a2bSJed Brown   c->free_ij = PETSC_TRUE;
5940ced3a2bSJed Brown   c->nonew   = 0;
59526fbe8dcSKarl Rupp 
59689d95c1aSJed Brown   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ;
5970ced3a2bSJed Brown 
5980ced3a2bSJed Brown   /* set MatInfo */
5990ced3a2bSJed Brown   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
6000ced3a2bSJed Brown   if (afill < 1.0) afill = 1.0;
6010ced3a2bSJed Brown   c->maxnz                     = ci[am];
6020ced3a2bSJed Brown   c->nz                        = ci[am];
6030ced3a2bSJed Brown   (*C)->info.mallocs           = ndouble;
6040ced3a2bSJed Brown   (*C)->info.fill_ratio_given  = fill;
6050ced3a2bSJed Brown   (*C)->info.fill_ratio_needed = afill;
6060ced3a2bSJed Brown 
6070ced3a2bSJed Brown #if defined(PETSC_USE_INFO)
6080ced3a2bSJed Brown   if (ci[am]) {
60957622a8eSBarry Smith     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %g needed %g.\n",ndouble,(double)fill,(double)afill);CHKERRQ(ierr);
61057622a8eSBarry Smith     ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr);
6110ced3a2bSJed Brown   } else {
6120ced3a2bSJed Brown     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
6130ced3a2bSJed Brown   }
6140ced3a2bSJed Brown #endif
6150ced3a2bSJed Brown   PetscFunctionReturn(0);
6160ced3a2bSJed Brown }
617e9e4536cSHong Zhang 
6188a07c6f1SJed Brown PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_BTHeap(Mat A,Mat B,PetscReal fill,Mat *C)
6198a07c6f1SJed Brown {
6208a07c6f1SJed Brown   PetscErrorCode     ierr;
6218a07c6f1SJed Brown   Mat_SeqAIJ         *a  = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
6228a07c6f1SJed Brown   const PetscInt     *ai = a->i,*bi=b->i,*aj=a->j,*bj=b->j;
6238a07c6f1SJed Brown   PetscInt           *ci,*cj,*bb;
6248a07c6f1SJed Brown   PetscInt           am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
6258a07c6f1SJed Brown   PetscReal          afill;
6268a07c6f1SJed Brown   PetscInt           i,j,col,ndouble = 0;
6270298fd71SBarry Smith   PetscFreeSpaceList free_space=NULL,current_space=NULL;
6288a07c6f1SJed Brown   PetscHeap          h;
6298a07c6f1SJed Brown   PetscBT            bt;
6308a07c6f1SJed Brown 
6318a07c6f1SJed Brown   PetscFunctionBegin;
632cd093f1dSJed Brown   /* Get ci and cj - using a heap for the sorted rows, but use BT so that each index is only added once */
6338a07c6f1SJed Brown   /*---------------------------------------------------------------------------------------------*/
6348a07c6f1SJed Brown   /* Allocate arrays for fill computation and free space for accumulating nonzero column */
635854ce69bSBarry Smith   ierr  = PetscMalloc1(am+2,&ci);CHKERRQ(ierr);
6368a07c6f1SJed Brown   ci[0] = 0;
6378a07c6f1SJed Brown 
6388a07c6f1SJed Brown   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
639f91af8c7SBarry Smith   ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ai[am],bi[bm])),&free_space);CHKERRQ(ierr);
6402205254eSKarl Rupp 
6418a07c6f1SJed Brown   current_space = free_space;
6428a07c6f1SJed Brown 
6438a07c6f1SJed Brown   ierr = PetscHeapCreate(a->rmax,&h);CHKERRQ(ierr);
644785e854fSJed Brown   ierr = PetscMalloc1(a->rmax,&bb);CHKERRQ(ierr);
6458a07c6f1SJed Brown   ierr = PetscBTCreate(bn,&bt);CHKERRQ(ierr);
6468a07c6f1SJed Brown 
6478a07c6f1SJed Brown   /* Determine ci and cj */
6488a07c6f1SJed Brown   for (i=0; i<am; i++) {
6498a07c6f1SJed 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 */
6508a07c6f1SJed Brown     const PetscInt *acol = aj + ai[i]; /* column indices of nonzero entries in this row */
6518a07c6f1SJed Brown     const PetscInt *fptr = current_space->array; /* Save beginning of the row so we can clear the BT later */
6528a07c6f1SJed Brown     ci[i+1] = ci[i];
6538a07c6f1SJed Brown     /* Populate the min heap */
6548a07c6f1SJed Brown     for (j=0; j<anzi; j++) {
6558a07c6f1SJed Brown       PetscInt brow = acol[j];
6568a07c6f1SJed Brown       for (bb[j] = bi[brow]; bb[j] < bi[brow+1]; bb[j]++) {
6578a07c6f1SJed Brown         PetscInt bcol = bj[bb[j]];
6588a07c6f1SJed Brown         if (!PetscBTLookupSet(bt,bcol)) { /* new entry */
6598a07c6f1SJed Brown           ierr = PetscHeapAdd(h,j,bcol);CHKERRQ(ierr);
6608a07c6f1SJed Brown           bb[j]++;
6618a07c6f1SJed Brown           break;
6628a07c6f1SJed Brown         }
6638a07c6f1SJed Brown       }
6648a07c6f1SJed Brown     }
6658a07c6f1SJed Brown     /* Pick off the min element, adding it to free space */
6668a07c6f1SJed Brown     ierr = PetscHeapPop(h,&j,&col);CHKERRQ(ierr);
6678a07c6f1SJed Brown     while (j >= 0) {
6688a07c6f1SJed Brown       if (current_space->local_remaining < 1) { /* double the size, but don't exceed 16 MiB */
6690298fd71SBarry Smith         fptr = NULL;                      /* need PetscBTMemzero */
670f91af8c7SBarry Smith         ierr = PetscFreeSpaceGet(PetscMin(PetscIntMultTruncate(2,current_space->total_array_size),16 << 20),&current_space);CHKERRQ(ierr);
6718a07c6f1SJed Brown         ndouble++;
6728a07c6f1SJed Brown       }
6738a07c6f1SJed Brown       *(current_space->array++) = col;
6748a07c6f1SJed Brown       current_space->local_used++;
6758a07c6f1SJed Brown       current_space->local_remaining--;
6768a07c6f1SJed Brown       ci[i+1]++;
6778a07c6f1SJed Brown 
6788a07c6f1SJed Brown       /* stash if anything else remains in this row of B */
6798a07c6f1SJed Brown       for (; bb[j] < bi[acol[j]+1]; bb[j]++) {
6808a07c6f1SJed Brown         PetscInt bcol = bj[bb[j]];
6818a07c6f1SJed Brown         if (!PetscBTLookupSet(bt,bcol)) { /* new entry */
6828a07c6f1SJed Brown           ierr = PetscHeapAdd(h,j,bcol);CHKERRQ(ierr);
6838a07c6f1SJed Brown           bb[j]++;
6848a07c6f1SJed Brown           break;
6858a07c6f1SJed Brown         }
6868a07c6f1SJed Brown       }
6878a07c6f1SJed Brown       ierr = PetscHeapPop(h,&j,&col);CHKERRQ(ierr);
6888a07c6f1SJed Brown     }
6898a07c6f1SJed Brown     if (fptr) {                 /* Clear the bits for this row */
6908a07c6f1SJed Brown       for (; fptr<current_space->array; fptr++) {ierr = PetscBTClear(bt,*fptr);CHKERRQ(ierr);}
6918a07c6f1SJed Brown     } else {                    /* We reallocated so we don't remember (easily) how to clear only the bits we changed */
6928a07c6f1SJed Brown       ierr = PetscBTMemzero(bn,bt);CHKERRQ(ierr);
6938a07c6f1SJed Brown     }
6948a07c6f1SJed Brown   }
6958a07c6f1SJed Brown   ierr = PetscFree(bb);CHKERRQ(ierr);
6968a07c6f1SJed Brown   ierr = PetscHeapDestroy(&h);CHKERRQ(ierr);
6978a07c6f1SJed Brown   ierr = PetscBTDestroy(&bt);CHKERRQ(ierr);
6988a07c6f1SJed Brown 
6998a07c6f1SJed Brown   /* Column indices are in the list of free space */
7008a07c6f1SJed Brown   /* Allocate space for cj, initialize cj, and */
7018a07c6f1SJed Brown   /* destroy list of free space and other temporary array(s) */
702785e854fSJed Brown   ierr = PetscMalloc1(ci[am],&cj);CHKERRQ(ierr);
7038a07c6f1SJed Brown   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
7048a07c6f1SJed Brown 
7058a07c6f1SJed Brown   /* put together the new symbolic matrix */
706ce94432eSBarry Smith   ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,NULL,C);CHKERRQ(ierr);
70733d57670SJed Brown   ierr = MatSetBlockSizesFromMats(*C,A,B);CHKERRQ(ierr);
708*02fe1965SBarry Smith   ierr = MatSetType(*C,((PetscObject)A)->type_name);CHKERRQ(ierr);
7098a07c6f1SJed Brown 
7108a07c6f1SJed Brown   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
7118a07c6f1SJed Brown   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
7128a07c6f1SJed Brown   c          = (Mat_SeqAIJ*)((*C)->data);
7138a07c6f1SJed Brown   c->free_a  = PETSC_TRUE;
7148a07c6f1SJed Brown   c->free_ij = PETSC_TRUE;
7158a07c6f1SJed Brown   c->nonew   = 0;
71626fbe8dcSKarl Rupp 
71789d95c1aSJed Brown   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ;
7188a07c6f1SJed Brown 
7198a07c6f1SJed Brown   /* set MatInfo */
7208a07c6f1SJed Brown   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
7218a07c6f1SJed Brown   if (afill < 1.0) afill = 1.0;
7228a07c6f1SJed Brown   c->maxnz                     = ci[am];
7238a07c6f1SJed Brown   c->nz                        = ci[am];
7248a07c6f1SJed Brown   (*C)->info.mallocs           = ndouble;
7258a07c6f1SJed Brown   (*C)->info.fill_ratio_given  = fill;
7268a07c6f1SJed Brown   (*C)->info.fill_ratio_needed = afill;
7278a07c6f1SJed Brown 
7288a07c6f1SJed Brown #if defined(PETSC_USE_INFO)
7298a07c6f1SJed Brown   if (ci[am]) {
73057622a8eSBarry Smith     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %g needed %g.\n",ndouble,(double)fill,(double)afill);CHKERRQ(ierr);
73157622a8eSBarry Smith     ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr);
7328a07c6f1SJed Brown   } else {
7338a07c6f1SJed Brown     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
7348a07c6f1SJed Brown   }
7358a07c6f1SJed Brown #endif
7368a07c6f1SJed Brown   PetscFunctionReturn(0);
7378a07c6f1SJed Brown }
7388a07c6f1SJed Brown 
739cd093f1dSJed Brown /* concatenate unique entries and then sort */
74058cf0668SJed Brown PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
741cd093f1dSJed Brown {
742cd093f1dSJed Brown   PetscErrorCode     ierr;
743cd093f1dSJed Brown   Mat_SeqAIJ         *a  = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
744cd093f1dSJed Brown   const PetscInt     *ai = a->i,*bi=b->i,*aj=a->j,*bj=b->j;
745cd093f1dSJed Brown   PetscInt           *ci,*cj;
746cd093f1dSJed Brown   PetscInt           am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
747cd093f1dSJed Brown   PetscReal          afill;
748cd093f1dSJed Brown   PetscInt           i,j,ndouble = 0;
749cd093f1dSJed Brown   PetscSegBuffer     seg,segrow;
750cd093f1dSJed Brown   char               *seen;
751cd093f1dSJed Brown 
752cd093f1dSJed Brown   PetscFunctionBegin;
753854ce69bSBarry Smith   ierr  = PetscMalloc1(am+1,&ci);CHKERRQ(ierr);
754cd093f1dSJed Brown   ci[0] = 0;
755cd093f1dSJed Brown 
756cd093f1dSJed Brown   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
757cd093f1dSJed Brown   ierr = PetscSegBufferCreate(sizeof(PetscInt),(PetscInt)(fill*(ai[am]+bi[bm])),&seg);CHKERRQ(ierr);
758cd093f1dSJed Brown   ierr = PetscSegBufferCreate(sizeof(PetscInt),100,&segrow);CHKERRQ(ierr);
759785e854fSJed Brown   ierr = PetscMalloc1(bn,&seen);CHKERRQ(ierr);
760cd093f1dSJed Brown   ierr = PetscMemzero(seen,bn*sizeof(char));CHKERRQ(ierr);
761cd093f1dSJed Brown 
762cd093f1dSJed Brown   /* Determine ci and cj */
763cd093f1dSJed Brown   for (i=0; i<am; i++) {
764cd093f1dSJed 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 */
765cd093f1dSJed Brown     const PetscInt *acol = aj + ai[i]; /* column indices of nonzero entries in this row */
766cd093f1dSJed Brown     PetscInt packlen = 0,*PETSC_RESTRICT crow;
767cd093f1dSJed Brown     /* Pack segrow */
768cd093f1dSJed Brown     for (j=0; j<anzi; j++) {
769cd093f1dSJed Brown       PetscInt brow = acol[j],bjstart = bi[brow],bjend = bi[brow+1],k;
770cd093f1dSJed Brown       for (k=bjstart; k<bjend; k++) {
771cd093f1dSJed Brown         PetscInt bcol = bj[k];
772cd093f1dSJed Brown         if (!seen[bcol]) { /* new entry */
773cd093f1dSJed Brown           PetscInt *PETSC_RESTRICT slot;
774cd093f1dSJed Brown           ierr = PetscSegBufferGetInts(segrow,1,&slot);CHKERRQ(ierr);
775cd093f1dSJed Brown           *slot = bcol;
776cd093f1dSJed Brown           seen[bcol] = 1;
777cd093f1dSJed Brown           packlen++;
778cd093f1dSJed Brown         }
779cd093f1dSJed Brown       }
780cd093f1dSJed Brown     }
781cd093f1dSJed Brown     ierr = PetscSegBufferGetInts(seg,packlen,&crow);CHKERRQ(ierr);
782cd093f1dSJed Brown     ierr = PetscSegBufferExtractTo(segrow,crow);CHKERRQ(ierr);
783cd093f1dSJed Brown     ierr = PetscSortInt(packlen,crow);CHKERRQ(ierr);
784cd093f1dSJed Brown     ci[i+1] = ci[i] + packlen;
785cd093f1dSJed Brown     for (j=0; j<packlen; j++) seen[crow[j]] = 0;
786cd093f1dSJed Brown   }
787cd093f1dSJed Brown   ierr = PetscSegBufferDestroy(&segrow);CHKERRQ(ierr);
788cd093f1dSJed Brown   ierr = PetscFree(seen);CHKERRQ(ierr);
789cd093f1dSJed Brown 
790cd093f1dSJed Brown   /* Column indices are in the segmented buffer */
791cd093f1dSJed Brown   ierr = PetscSegBufferExtractAlloc(seg,&cj);CHKERRQ(ierr);
792cd093f1dSJed Brown   ierr = PetscSegBufferDestroy(&seg);CHKERRQ(ierr);
793cd093f1dSJed Brown 
794cd093f1dSJed Brown   /* put together the new symbolic matrix */
795cd093f1dSJed Brown   ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,NULL,C);CHKERRQ(ierr);
79633d57670SJed Brown   ierr = MatSetBlockSizesFromMats(*C,A,B);CHKERRQ(ierr);
797*02fe1965SBarry Smith   ierr = MatSetType(*C,((PetscObject)A)->type_name);CHKERRQ(ierr);
798cd093f1dSJed Brown 
799cd093f1dSJed Brown   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
800cd093f1dSJed Brown   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
801cd093f1dSJed Brown   c          = (Mat_SeqAIJ*)((*C)->data);
802cd093f1dSJed Brown   c->free_a  = PETSC_TRUE;
803cd093f1dSJed Brown   c->free_ij = PETSC_TRUE;
804cd093f1dSJed Brown   c->nonew   = 0;
805cd093f1dSJed Brown 
806cd093f1dSJed Brown   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ;
807cd093f1dSJed Brown 
808cd093f1dSJed Brown   /* set MatInfo */
809cd093f1dSJed Brown   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
810cd093f1dSJed Brown   if (afill < 1.0) afill = 1.0;
811cd093f1dSJed Brown   c->maxnz                     = ci[am];
812cd093f1dSJed Brown   c->nz                        = ci[am];
813cd093f1dSJed Brown   (*C)->info.mallocs           = ndouble;
814cd093f1dSJed Brown   (*C)->info.fill_ratio_given  = fill;
815cd093f1dSJed Brown   (*C)->info.fill_ratio_needed = afill;
816cd093f1dSJed Brown 
817cd093f1dSJed Brown #if defined(PETSC_USE_INFO)
818cd093f1dSJed Brown   if (ci[am]) {
81957622a8eSBarry Smith     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %g needed %g.\n",ndouble,(double)fill,(double)afill);CHKERRQ(ierr);
82057622a8eSBarry Smith     ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr);
821cd093f1dSJed Brown   } else {
822cd093f1dSJed Brown     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
823cd093f1dSJed Brown   }
824cd093f1dSJed Brown #endif
825cd093f1dSJed Brown   PetscFunctionReturn(0);
826cd093f1dSJed Brown }
827cd093f1dSJed Brown 
828d2853540SHong Zhang /* This routine is not used. Should be removed! */
8296fc122caSHong Zhang PetscErrorCode MatMatTransposeMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
8305df89d91SHong Zhang {
831bc011b1eSHong Zhang   PetscErrorCode ierr;
832bc011b1eSHong Zhang 
833bc011b1eSHong Zhang   PetscFunctionBegin;
834bc011b1eSHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
8353ff4c91cSHong Zhang     ierr = PetscLogEventBegin(MAT_MatTransposeMultSymbolic,A,B,0,0);CHKERRQ(ierr);
8366fc122caSHong Zhang     ierr = MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr);
8373ff4c91cSHong Zhang     ierr = PetscLogEventEnd(MAT_MatTransposeMultSymbolic,A,B,0,0);CHKERRQ(ierr);
838bc011b1eSHong Zhang   }
8393ff4c91cSHong Zhang   ierr = PetscLogEventBegin(MAT_MatTransposeMultNumeric,A,B,0,0);CHKERRQ(ierr);
8406fc122caSHong Zhang   ierr = MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr);
8413ff4c91cSHong Zhang   ierr = PetscLogEventEnd(MAT_MatTransposeMultNumeric,A,B,0,0);CHKERRQ(ierr);
842bc011b1eSHong Zhang   PetscFunctionReturn(0);
843bc011b1eSHong Zhang }
844bc011b1eSHong Zhang 
8452128a86cSHong Zhang PetscErrorCode MatDestroy_SeqAIJ_MatMatMultTrans(Mat A)
8462128a86cSHong Zhang {
8472128a86cSHong Zhang   PetscErrorCode      ierr;
8484c7df5ccSHong Zhang   Mat_SeqAIJ          *a=(Mat_SeqAIJ*)A->data;
84940192850SHong Zhang   Mat_MatMatTransMult *abt=a->abt;
8502128a86cSHong Zhang 
8512128a86cSHong Zhang   PetscFunctionBegin;
85240192850SHong Zhang   ierr = (abt->destroy)(A);CHKERRQ(ierr);
85340192850SHong Zhang   ierr = MatTransposeColoringDestroy(&abt->matcoloring);CHKERRQ(ierr);
85440192850SHong Zhang   ierr = MatDestroy(&abt->Bt_den);CHKERRQ(ierr);
85540192850SHong Zhang   ierr = MatDestroy(&abt->ABt_den);CHKERRQ(ierr);
85640192850SHong Zhang   ierr = PetscFree(abt);CHKERRQ(ierr);
8572128a86cSHong Zhang   PetscFunctionReturn(0);
8582128a86cSHong Zhang }
8592128a86cSHong Zhang 
8606fc122caSHong Zhang PetscErrorCode MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
861bc011b1eSHong Zhang {
8625df89d91SHong Zhang   PetscErrorCode      ierr;
86381d82fe4SHong Zhang   Mat                 Bt;
86481d82fe4SHong Zhang   PetscInt            *bti,*btj;
86540192850SHong Zhang   Mat_MatMatTransMult *abt;
8664c7df5ccSHong Zhang   Mat_SeqAIJ          *c;
867d2853540SHong Zhang 
86881d82fe4SHong Zhang   PetscFunctionBegin;
86981d82fe4SHong Zhang   /* create symbolic Bt */
87081d82fe4SHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(B,&bti,&btj);CHKERRQ(ierr);
8710298fd71SBarry Smith   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,B->cmap->n,B->rmap->n,bti,btj,NULL,&Bt);CHKERRQ(ierr);
87233d57670SJed Brown   ierr = MatSetBlockSizes(Bt,PetscAbs(A->cmap->bs),PetscAbs(B->cmap->bs));CHKERRQ(ierr);
873*02fe1965SBarry Smith   ierr = MatSetType(*C,((PetscObject)A)->type_name);CHKERRQ(ierr);
87481d82fe4SHong Zhang 
87581d82fe4SHong Zhang   /* get symbolic C=A*Bt */
87681d82fe4SHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,Bt,fill,C);CHKERRQ(ierr);
87781d82fe4SHong Zhang 
8782128a86cSHong Zhang   /* create a supporting struct for reuse intermidiate dense matrices with matcoloring */
879b00a9115SJed Brown   ierr   = PetscNew(&abt);CHKERRQ(ierr);
8804c7df5ccSHong Zhang   c      = (Mat_SeqAIJ*)(*C)->data;
88140192850SHong Zhang   c->abt = abt;
8822128a86cSHong Zhang 
88340192850SHong Zhang   abt->usecoloring = PETSC_FALSE;
88440192850SHong Zhang   abt->destroy     = (*C)->ops->destroy;
8852128a86cSHong Zhang   (*C)->ops->destroy     = MatDestroy_SeqAIJ_MatMatMultTrans;
8862128a86cSHong Zhang 
887c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)A)->options,NULL,"-matmattransmult_color",&abt->usecoloring,NULL);CHKERRQ(ierr);
88840192850SHong Zhang   if (abt->usecoloring) {
889b9af6bddSHong Zhang     /* Create MatTransposeColoring from symbolic C=A*B^T */
890b9af6bddSHong Zhang     MatTransposeColoring matcoloring;
891335efc43SPeter Brune     MatColoring          coloring;
8922128a86cSHong Zhang     ISColoring           iscoloring;
8932128a86cSHong Zhang     Mat                  Bt_dense,C_dense;
8944d478ae7SHong Zhang     Mat_SeqAIJ           *c=(Mat_SeqAIJ*)(*C)->data;
8954d478ae7SHong Zhang     /* inode causes memory problem, don't know why */
8964d478ae7SHong Zhang     if (c->inode.use) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MAT_USE_INODES is not supported. Use '-mat_no_inode'");
897e8354b3eSHong Zhang 
898335efc43SPeter Brune     ierr = MatColoringCreate(*C,&coloring);CHKERRQ(ierr);
899335efc43SPeter Brune     ierr = MatColoringSetDistance(coloring,2);CHKERRQ(ierr);
900335efc43SPeter Brune     ierr = MatColoringSetType(coloring,MATCOLORINGSL);CHKERRQ(ierr);
901335efc43SPeter Brune     ierr = MatColoringSetFromOptions(coloring);CHKERRQ(ierr);
902335efc43SPeter Brune     ierr = MatColoringApply(coloring,&iscoloring);CHKERRQ(ierr);
903335efc43SPeter Brune     ierr = MatColoringDestroy(&coloring);CHKERRQ(ierr);
904b9af6bddSHong Zhang     ierr = MatTransposeColoringCreate(*C,iscoloring,&matcoloring);CHKERRQ(ierr);
9052205254eSKarl Rupp 
90640192850SHong Zhang     abt->matcoloring = matcoloring;
9072205254eSKarl Rupp 
9082128a86cSHong Zhang     ierr = ISColoringDestroy(&iscoloring);CHKERRQ(ierr);
9092128a86cSHong Zhang 
9102128a86cSHong Zhang     /* Create Bt_dense and C_dense = A*Bt_dense */
9112128a86cSHong Zhang     ierr = MatCreate(PETSC_COMM_SELF,&Bt_dense);CHKERRQ(ierr);
9122128a86cSHong Zhang     ierr = MatSetSizes(Bt_dense,A->cmap->n,matcoloring->ncolors,A->cmap->n,matcoloring->ncolors);CHKERRQ(ierr);
9132128a86cSHong Zhang     ierr = MatSetType(Bt_dense,MATSEQDENSE);CHKERRQ(ierr);
9140298fd71SBarry Smith     ierr = MatSeqDenseSetPreallocation(Bt_dense,NULL);CHKERRQ(ierr);
9152205254eSKarl Rupp 
9162128a86cSHong Zhang     Bt_dense->assembled = PETSC_TRUE;
91740192850SHong Zhang     abt->Bt_den   = Bt_dense;
9182128a86cSHong Zhang 
9192128a86cSHong Zhang     ierr = MatCreate(PETSC_COMM_SELF,&C_dense);CHKERRQ(ierr);
9202128a86cSHong Zhang     ierr = MatSetSizes(C_dense,A->rmap->n,matcoloring->ncolors,A->rmap->n,matcoloring->ncolors);CHKERRQ(ierr);
9212128a86cSHong Zhang     ierr = MatSetType(C_dense,MATSEQDENSE);CHKERRQ(ierr);
9220298fd71SBarry Smith     ierr = MatSeqDenseSetPreallocation(C_dense,NULL);CHKERRQ(ierr);
9232205254eSKarl Rupp 
9242128a86cSHong Zhang     Bt_dense->assembled = PETSC_TRUE;
92540192850SHong Zhang     abt->ABt_den  = C_dense;
926f94ccd6cSHong Zhang 
927f94ccd6cSHong Zhang #if defined(PETSC_USE_INFO)
9281ce71dffSSatish Balay     {
929f94ccd6cSHong Zhang       Mat_SeqAIJ *c = (Mat_SeqAIJ*)(*C)->data;
930c40ebe3bSHong 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);
9311ce71dffSSatish Balay     }
932f94ccd6cSHong Zhang #endif
9332128a86cSHong Zhang   }
934e99be685SHong Zhang   /* clean up */
935e99be685SHong Zhang   ierr = MatDestroy(&Bt);CHKERRQ(ierr);
936e99be685SHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(B,&bti,&btj);CHKERRQ(ierr);
9375df89d91SHong Zhang   PetscFunctionReturn(0);
9385df89d91SHong Zhang }
9395df89d91SHong Zhang 
9406fc122caSHong Zhang PetscErrorCode MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
9415df89d91SHong Zhang {
9425df89d91SHong Zhang   PetscErrorCode      ierr;
9435df89d91SHong Zhang   Mat_SeqAIJ          *a   =(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data;
944e2cac8caSJed Brown   PetscInt            *ai  =a->i,*aj=a->j,*bi=b->i,*bj=b->j,anzi,bnzj,nexta,nextb,*acol,*bcol,brow;
94581d82fe4SHong Zhang   PetscInt            cm   =C->rmap->n,*ci=c->i,*cj=c->j,i,j,cnzi,*ccol;
9465df89d91SHong Zhang   PetscLogDouble      flops=0.0;
947aa1aec99SHong Zhang   MatScalar           *aa  =a->a,*aval,*ba=b->a,*bval,*ca,*cval;
94840192850SHong Zhang   Mat_MatMatTransMult *abt = c->abt;
9495df89d91SHong Zhang 
9505df89d91SHong Zhang   PetscFunctionBegin;
95158ed3ceaSHong Zhang   /* clear old values in C */
95258ed3ceaSHong Zhang   if (!c->a) {
953854ce69bSBarry Smith     ierr      = PetscMalloc1(ci[cm]+1,&ca);CHKERRQ(ierr);
95458ed3ceaSHong Zhang     c->a      = ca;
95558ed3ceaSHong Zhang     c->free_a = PETSC_TRUE;
95658ed3ceaSHong Zhang   } else {
95758ed3ceaSHong Zhang     ca =  c->a;
95858ed3ceaSHong Zhang   }
95958ed3ceaSHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
96058ed3ceaSHong Zhang 
96140192850SHong Zhang   if (abt->usecoloring) {
96240192850SHong Zhang     MatTransposeColoring matcoloring = abt->matcoloring;
96340192850SHong Zhang     Mat                  Bt_dense,C_dense = abt->ABt_den;
964c8db22f6SHong Zhang 
965b9af6bddSHong Zhang     /* Get Bt_dense by Apply MatTransposeColoring to B */
96640192850SHong Zhang     Bt_dense = abt->Bt_den;
967b9af6bddSHong Zhang     ierr = MatTransColoringApplySpToDen(matcoloring,B,Bt_dense);CHKERRQ(ierr);
968c8db22f6SHong Zhang 
969c8db22f6SHong Zhang     /* C_dense = A*Bt_dense */
9702128a86cSHong Zhang     ierr = MatMatMultNumeric_SeqAIJ_SeqDense(A,Bt_dense,C_dense);CHKERRQ(ierr);
971c8db22f6SHong Zhang 
9722128a86cSHong Zhang     /* Recover C from C_dense */
973b9af6bddSHong Zhang     ierr = MatTransColoringApplyDenToSp(matcoloring,C_dense,C);CHKERRQ(ierr);
974c8db22f6SHong Zhang     PetscFunctionReturn(0);
975c8db22f6SHong Zhang   }
976c8db22f6SHong Zhang 
9771fa1209cSHong Zhang   for (i=0; i<cm; i++) {
97881d82fe4SHong Zhang     anzi = ai[i+1] - ai[i];
97981d82fe4SHong Zhang     acol = aj + ai[i];
9806973769bSHong Zhang     aval = aa + ai[i];
9811fa1209cSHong Zhang     cnzi = ci[i+1] - ci[i];
9821fa1209cSHong Zhang     ccol = cj + ci[i];
9836973769bSHong Zhang     cval = ca + ci[i];
9841fa1209cSHong Zhang     for (j=0; j<cnzi; j++) {
98581d82fe4SHong Zhang       brow = ccol[j];
98681d82fe4SHong Zhang       bnzj = bi[brow+1] - bi[brow];
98781d82fe4SHong Zhang       bcol = bj + bi[brow];
9886973769bSHong Zhang       bval = ba + bi[brow];
9896973769bSHong Zhang 
99081d82fe4SHong Zhang       /* perform sparse inner-product c(i,j)=A[i,:]*B[j,:]^T */
99181d82fe4SHong Zhang       nexta = 0; nextb = 0;
99281d82fe4SHong Zhang       while (nexta<anzi && nextb<bnzj) {
9937b6d5e96SMark Adams         while (nexta < anzi && acol[nexta] < bcol[nextb]) nexta++;
99481d82fe4SHong Zhang         if (nexta == anzi) break;
9957b6d5e96SMark Adams         while (nextb < bnzj && acol[nexta] > bcol[nextb]) nextb++;
99681d82fe4SHong Zhang         if (nextb == bnzj) break;
99781d82fe4SHong Zhang         if (acol[nexta] == bcol[nextb]) {
9986973769bSHong Zhang           cval[j] += aval[nexta]*bval[nextb];
99981d82fe4SHong Zhang           nexta++; nextb++;
100081d82fe4SHong Zhang           flops += 2;
10011fa1209cSHong Zhang         }
10021fa1209cSHong Zhang       }
100381d82fe4SHong Zhang     }
100481d82fe4SHong Zhang   }
100581d82fe4SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
100681d82fe4SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
100781d82fe4SHong Zhang   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
10085df89d91SHong Zhang   PetscFunctionReturn(0);
10095df89d91SHong Zhang }
10105df89d91SHong Zhang 
10116d373c3eSHong Zhang PetscErrorCode MatDestroy_SeqAIJ_MatTransMatMult(Mat A)
10126d373c3eSHong Zhang {
10136d373c3eSHong Zhang   PetscErrorCode      ierr;
10146d373c3eSHong Zhang   Mat_SeqAIJ          *a = (Mat_SeqAIJ*)A->data;
10156d373c3eSHong Zhang   Mat_MatTransMatMult *atb = a->atb;
10166d373c3eSHong Zhang 
10176d373c3eSHong Zhang   PetscFunctionBegin;
10186d373c3eSHong Zhang   ierr = MatDestroy(&atb->At);CHKERRQ(ierr);
10196d373c3eSHong Zhang   ierr = (atb->destroy)(A);CHKERRQ(ierr);
10206d373c3eSHong Zhang   ierr = PetscFree(atb);CHKERRQ(ierr);
10216d373c3eSHong Zhang   PetscFunctionReturn(0);
10226d373c3eSHong Zhang }
10236d373c3eSHong Zhang 
10240adebc6cSBarry Smith PetscErrorCode MatTransposeMatMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
10250adebc6cSBarry Smith {
10265df89d91SHong Zhang   PetscErrorCode      ierr;
10276d373c3eSHong Zhang   const char          *algTypes[2] = {"matmatmult","outerproduct"};
10286d373c3eSHong Zhang   PetscInt            alg=0; /* set default algorithm */
10296d373c3eSHong Zhang   Mat                 At;
10306d373c3eSHong Zhang   Mat_MatTransMatMult *atb;
10316d373c3eSHong Zhang   Mat_SeqAIJ          *c;
10325df89d91SHong Zhang 
10335df89d91SHong Zhang   PetscFunctionBegin;
10345df89d91SHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
10356d373c3eSHong Zhang     ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr);
10366d373c3eSHong Zhang     PetscOptionsObject->alreadyprinted = PETSC_FALSE; /* a hack to ensure the option shows in '-help' */
10376d373c3eSHong Zhang     ierr = PetscOptionsEList("-mattransposematmult_via","Algorithmic approach","MatTransposeMatMult",algTypes,2,algTypes[0],&alg,NULL);CHKERRQ(ierr);
10386d373c3eSHong Zhang     ierr = PetscOptionsEnd();CHKERRQ(ierr);
10396d373c3eSHong Zhang 
10406d373c3eSHong Zhang     switch (alg) {
10416d373c3eSHong Zhang     case 1:
104275648e8dSHong Zhang       ierr = MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr);
10436d373c3eSHong Zhang       break;
10446d373c3eSHong Zhang     default:
10456d373c3eSHong Zhang       ierr = PetscNew(&atb);CHKERRQ(ierr);
10466d373c3eSHong Zhang       ierr = MatTranspose_SeqAIJ(A,MAT_INITIAL_MATRIX,&At);CHKERRQ(ierr);
10476d373c3eSHong Zhang       ierr = MatMatMult_SeqAIJ_SeqAIJ(At,B,MAT_INITIAL_MATRIX,fill,C);CHKERRQ(ierr);
10486d373c3eSHong Zhang 
1049618cf492SHong Zhang       c                  = (Mat_SeqAIJ*)(*C)->data;
10506d373c3eSHong Zhang       c->atb             = atb;
10516d373c3eSHong Zhang       atb->At            = At;
10526d373c3eSHong Zhang       atb->destroy       = (*C)->ops->destroy;
10536d373c3eSHong Zhang       (*C)->ops->destroy = MatDestroy_SeqAIJ_MatTransMatMult;
10546d373c3eSHong Zhang 
10556d373c3eSHong Zhang       break;
10565df89d91SHong Zhang     }
10576d373c3eSHong Zhang   }
10586d373c3eSHong Zhang   if (alg) {
10596d373c3eSHong Zhang     ierr = (*(*C)->ops->mattransposemultnumeric)(A,B,*C);CHKERRQ(ierr);
10606d373c3eSHong Zhang   } else if (!alg && scall == MAT_REUSE_MATRIX) {
10616d373c3eSHong Zhang     c   = (Mat_SeqAIJ*)(*C)->data;
10626d373c3eSHong Zhang     atb = c->atb;
10636d373c3eSHong Zhang     At  = atb->At;
10646d373c3eSHong Zhang     ierr = MatTranspose_SeqAIJ(A,MAT_REUSE_MATRIX,&At);CHKERRQ(ierr);
10656d373c3eSHong Zhang     ierr = MatMatMult_SeqAIJ_SeqAIJ(At,B,MAT_REUSE_MATRIX,fill,C);CHKERRQ(ierr);
10666d373c3eSHong Zhang   }
10675df89d91SHong Zhang   PetscFunctionReturn(0);
10685df89d91SHong Zhang }
10695df89d91SHong Zhang 
107075648e8dSHong Zhang PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
10715df89d91SHong Zhang {
1072bc011b1eSHong Zhang   PetscErrorCode ierr;
1073bc011b1eSHong Zhang   Mat            At;
107438baddfdSBarry Smith   PetscInt       *ati,*atj;
1075bc011b1eSHong Zhang 
1076bc011b1eSHong Zhang   PetscFunctionBegin;
1077bc011b1eSHong Zhang   /* create symbolic At */
1078bc011b1eSHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr);
10790298fd71SBarry Smith   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,A->cmap->n,A->rmap->n,ati,atj,NULL,&At);CHKERRQ(ierr);
108033d57670SJed Brown   ierr = MatSetBlockSizes(At,PetscAbs(A->cmap->bs),PetscAbs(B->cmap->bs));CHKERRQ(ierr);
1081*02fe1965SBarry Smith   ierr = MatSetType(*C,((PetscObject)A)->type_name);CHKERRQ(ierr);
1082bc011b1eSHong Zhang 
1083bc011b1eSHong Zhang   /* get symbolic C=At*B */
1084bc011b1eSHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(At,B,fill,C);CHKERRQ(ierr);
1085bc011b1eSHong Zhang 
1086bc011b1eSHong Zhang   /* clean up */
10876bf464f9SBarry Smith   ierr = MatDestroy(&At);CHKERRQ(ierr);
1088bc011b1eSHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr);
10896d373c3eSHong Zhang 
10906d373c3eSHong Zhang   (*C)->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ;
1091bc011b1eSHong Zhang   PetscFunctionReturn(0);
1092bc011b1eSHong Zhang }
1093bc011b1eSHong Zhang 
109475648e8dSHong Zhang PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
1095bc011b1eSHong Zhang {
1096bc011b1eSHong Zhang   PetscErrorCode ierr;
10970fbc74f4SHong Zhang   Mat_SeqAIJ     *a   =(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data;
1098d0f46423SBarry Smith   PetscInt       am   =A->rmap->n,anzi,*ai=a->i,*aj=a->j,*bi=b->i,*bj,bnzi,nextb;
1099d0f46423SBarry Smith   PetscInt       cm   =C->rmap->n,*ci=c->i,*cj=c->j,crow,*cjj,i,j,k;
1100d13dce4bSSatish Balay   PetscLogDouble flops=0.0;
1101aa1aec99SHong Zhang   MatScalar      *aa  =a->a,*ba,*ca,*caj;
1102bc011b1eSHong Zhang 
1103bc011b1eSHong Zhang   PetscFunctionBegin;
1104aa1aec99SHong Zhang   if (!c->a) {
1105854ce69bSBarry Smith     ierr = PetscMalloc1(ci[cm]+1,&ca);CHKERRQ(ierr);
11062205254eSKarl Rupp 
1107aa1aec99SHong Zhang     c->a      = ca;
1108aa1aec99SHong Zhang     c->free_a = PETSC_TRUE;
1109aa1aec99SHong Zhang   } else {
1110aa1aec99SHong Zhang     ca = c->a;
1111aa1aec99SHong Zhang   }
1112bc011b1eSHong Zhang   /* clear old values in C */
1113bc011b1eSHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
1114bc011b1eSHong Zhang 
1115bc011b1eSHong Zhang   /* compute A^T*B using outer product (A^T)[:,i]*B[i,:] */
1116bc011b1eSHong Zhang   for (i=0; i<am; i++) {
1117bc011b1eSHong Zhang     bj   = b->j + bi[i];
1118bc011b1eSHong Zhang     ba   = b->a + bi[i];
1119bc011b1eSHong Zhang     bnzi = bi[i+1] - bi[i];
1120bc011b1eSHong Zhang     anzi = ai[i+1] - ai[i];
1121bc011b1eSHong Zhang     for (j=0; j<anzi; j++) {
1122bc011b1eSHong Zhang       nextb = 0;
11230fbc74f4SHong Zhang       crow  = *aj++;
1124bc011b1eSHong Zhang       cjj   = cj + ci[crow];
1125bc011b1eSHong Zhang       caj   = ca + ci[crow];
1126bc011b1eSHong Zhang       /* perform sparse axpy operation.  Note cjj includes bj. */
1127bc011b1eSHong Zhang       for (k=0; nextb<bnzi; k++) {
11280fbc74f4SHong Zhang         if (cjj[k] == *(bj+nextb)) { /* ccol == bcol */
11290fbc74f4SHong Zhang           caj[k] += (*aa)*(*(ba+nextb));
1130bc011b1eSHong Zhang           nextb++;
1131bc011b1eSHong Zhang         }
1132bc011b1eSHong Zhang       }
1133bc011b1eSHong Zhang       flops += 2*bnzi;
11340fbc74f4SHong Zhang       aa++;
1135bc011b1eSHong Zhang     }
1136bc011b1eSHong Zhang   }
1137bc011b1eSHong Zhang 
1138bc011b1eSHong Zhang   /* Assemble the final matrix and clean up */
1139bc011b1eSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1140bc011b1eSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1141bc011b1eSHong Zhang   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
1142bc011b1eSHong Zhang   PetscFunctionReturn(0);
1143bc011b1eSHong Zhang }
11449123193aSHong Zhang 
1145150d2497SBarry Smith PETSC_INTERN PetscErrorCode MatMatMult_SeqAIJ_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
11469123193aSHong Zhang {
11479123193aSHong Zhang   PetscErrorCode ierr;
11489123193aSHong Zhang 
11499123193aSHong Zhang   PetscFunctionBegin;
11509123193aSHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
11513ff4c91cSHong Zhang     ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
11529123193aSHong Zhang     ierr = MatMatMultSymbolic_SeqAIJ_SeqDense(A,B,fill,C);CHKERRQ(ierr);
11533ff4c91cSHong Zhang     ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
11549123193aSHong Zhang   }
11553ff4c91cSHong Zhang   ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
11569123193aSHong Zhang   ierr = MatMatMultNumeric_SeqAIJ_SeqDense(A,B,*C);CHKERRQ(ierr);
11574614e006SHong Zhang   ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
11589123193aSHong Zhang   PetscFunctionReturn(0);
11599123193aSHong Zhang }
1160edd81eecSMatthew Knepley 
11619123193aSHong Zhang PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C)
11629123193aSHong Zhang {
11639123193aSHong Zhang   PetscErrorCode ierr;
11649123193aSHong Zhang 
11659123193aSHong Zhang   PetscFunctionBegin;
11665a586d82SBarry Smith   ierr = MatMatMultSymbolic_SeqDense_SeqDense(A,B,0.0,C);CHKERRQ(ierr);
11672205254eSKarl Rupp 
1168d73949e8SHong Zhang   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqDense;
11699123193aSHong Zhang   PetscFunctionReturn(0);
11709123193aSHong Zhang }
11719123193aSHong Zhang 
11729123193aSHong Zhang PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqDense(Mat A,Mat B,Mat C)
11739123193aSHong Zhang {
1174f32d5d43SBarry Smith   Mat_SeqAIJ        *a=(Mat_SeqAIJ*)A->data;
1175612438f1SStefano Zampini   Mat_SeqDense      *bd = (Mat_SeqDense*)B->data;
11769123193aSHong Zhang   PetscErrorCode    ierr;
1177daf57211SHong Zhang   PetscScalar       *c,*b,r1,r2,r3,r4,*c1,*c2,*c3,*c4,aatmp;
1178daf57211SHong Zhang   const PetscScalar *aa,*b1,*b2,*b3,*b4;
1179daf57211SHong Zhang   const PetscInt    *aj;
1180612438f1SStefano Zampini   PetscInt          cm=C->rmap->n,cn=B->cmap->n,bm=bd->lda,am=A->rmap->n;
1181daf57211SHong Zhang   PetscInt          am4=4*am,bm4=4*bm,col,i,j,n,ajtmp;
11829123193aSHong Zhang 
11839123193aSHong Zhang   PetscFunctionBegin;
1184f32d5d43SBarry Smith   if (!cm || !cn) PetscFunctionReturn(0);
1185612438f1SStefano 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);
1186e32f2f54SBarry 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);
1187e32f2f54SBarry 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);
1188612438f1SStefano Zampini   b = bd->v;
11898c778c55SBarry Smith   ierr = MatDenseGetArray(C,&c);CHKERRQ(ierr);
1190f32d5d43SBarry Smith   b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm;
1191730858b9SHong Zhang   c1 = c; c2 = c1 + am; c3 = c2 + am; c4 = c3 + am;
1192f32d5d43SBarry Smith   for (col=0; col<cn-4; col += 4) {  /* over columns of C */
1193f32d5d43SBarry Smith     for (i=0; i<am; i++) {        /* over rows of C in those columns */
1194f32d5d43SBarry Smith       r1 = r2 = r3 = r4 = 0.0;
1195f32d5d43SBarry Smith       n  = a->i[i+1] - a->i[i];
1196f32d5d43SBarry Smith       aj = a->j + a->i[i];
1197f32d5d43SBarry Smith       aa = a->a + a->i[i];
1198f32d5d43SBarry Smith       for (j=0; j<n; j++) {
1199730858b9SHong Zhang         aatmp = aa[j]; ajtmp = aj[j];
1200730858b9SHong Zhang         r1 += aatmp*b1[ajtmp];
1201730858b9SHong Zhang         r2 += aatmp*b2[ajtmp];
1202730858b9SHong Zhang         r3 += aatmp*b3[ajtmp];
1203730858b9SHong Zhang         r4 += aatmp*b4[ajtmp];
12049123193aSHong Zhang       }
1205730858b9SHong Zhang       c1[i] = r1;
1206730858b9SHong Zhang       c2[i] = r2;
1207730858b9SHong Zhang       c3[i] = r3;
1208730858b9SHong Zhang       c4[i] = r4;
1209f32d5d43SBarry Smith     }
1210730858b9SHong Zhang     b1 += bm4; b2 += bm4; b3 += bm4; b4 += bm4;
1211730858b9SHong Zhang     c1 += am4; c2 += am4; c3 += am4; c4 += am4;
1212f32d5d43SBarry Smith   }
1213f32d5d43SBarry Smith   for (; col<cn; col++) {   /* over extra columns of C */
1214f32d5d43SBarry Smith     for (i=0; i<am; i++) {  /* over rows of C in those columns */
1215f32d5d43SBarry Smith       r1 = 0.0;
1216f32d5d43SBarry Smith       n  = a->i[i+1] - a->i[i];
1217f32d5d43SBarry Smith       aj = a->j + a->i[i];
1218f32d5d43SBarry Smith       aa = a->a + a->i[i];
1219f32d5d43SBarry Smith       for (j=0; j<n; j++) {
1220730858b9SHong Zhang         r1 += aa[j]*b1[aj[j]];
1221f32d5d43SBarry Smith       }
1222730858b9SHong Zhang       c1[i] = r1;
1223f32d5d43SBarry Smith     }
1224f32d5d43SBarry Smith     b1 += bm;
1225730858b9SHong Zhang     c1 += am;
1226f32d5d43SBarry Smith   }
1227dc0b31edSSatish Balay   ierr = PetscLogFlops(cn*(2.0*a->nz));CHKERRQ(ierr);
12288c778c55SBarry Smith   ierr = MatDenseRestoreArray(C,&c);CHKERRQ(ierr);
1229f32d5d43SBarry Smith   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1230f32d5d43SBarry Smith   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1231f32d5d43SBarry Smith   PetscFunctionReturn(0);
1232f32d5d43SBarry Smith }
1233f32d5d43SBarry Smith 
1234f32d5d43SBarry Smith /*
12354324174eSBarry Smith    Note very similar to MatMult_SeqAIJ(), should generate both codes from same base
1236f32d5d43SBarry Smith */
1237f32d5d43SBarry Smith PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat A,Mat B,Mat C)
1238f32d5d43SBarry Smith {
1239f32d5d43SBarry Smith   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
124090f5ea3eSStefano Zampini   Mat_SeqDense   *bd = (Mat_SeqDense*)B->data;
1241f32d5d43SBarry Smith   PetscErrorCode ierr;
1242dd6ea824SBarry Smith   PetscScalar    *b,*c,r1,r2,r3,r4,*b1,*b2,*b3,*b4;
1243dd6ea824SBarry Smith   MatScalar      *aa;
124490f5ea3eSStefano Zampini   PetscInt       cm  = C->rmap->n, cn=B->cmap->n, bm=bd->lda, col, i,j,n,*aj, am = A->rmap->n,*ii,arm;
12454324174eSBarry Smith   PetscInt       am2 = 2*am, am3 = 3*am,  bm4 = 4*bm,colam,*ridx;
1246f32d5d43SBarry Smith 
1247f32d5d43SBarry Smith   PetscFunctionBegin;
1248f32d5d43SBarry Smith   if (!cm || !cn) PetscFunctionReturn(0);
124990f5ea3eSStefano Zampini   b = bd->v;
12508c778c55SBarry Smith   ierr = MatDenseGetArray(C,&c);CHKERRQ(ierr);
1251f32d5d43SBarry Smith   b1   = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm;
12524324174eSBarry Smith 
12534324174eSBarry Smith   if (a->compressedrow.use) { /* use compressed row format */
12544324174eSBarry Smith     for (col=0; col<cn-4; col += 4) {  /* over columns of C */
12554324174eSBarry Smith       colam = col*am;
12564324174eSBarry Smith       arm   = a->compressedrow.nrows;
12574324174eSBarry Smith       ii    = a->compressedrow.i;
12584324174eSBarry Smith       ridx  = a->compressedrow.rindex;
12594324174eSBarry Smith       for (i=0; i<arm; i++) {        /* over rows of C in those columns */
12604324174eSBarry Smith         r1 = r2 = r3 = r4 = 0.0;
12614324174eSBarry Smith         n  = ii[i+1] - ii[i];
12624324174eSBarry Smith         aj = a->j + ii[i];
12634324174eSBarry Smith         aa = a->a + ii[i];
12644324174eSBarry Smith         for (j=0; j<n; j++) {
12654324174eSBarry Smith           r1 += (*aa)*b1[*aj];
12664324174eSBarry Smith           r2 += (*aa)*b2[*aj];
12674324174eSBarry Smith           r3 += (*aa)*b3[*aj];
12684324174eSBarry Smith           r4 += (*aa++)*b4[*aj++];
12694324174eSBarry Smith         }
12704324174eSBarry Smith         c[colam       + ridx[i]] += r1;
12714324174eSBarry Smith         c[colam + am  + ridx[i]] += r2;
12724324174eSBarry Smith         c[colam + am2 + ridx[i]] += r3;
12734324174eSBarry Smith         c[colam + am3 + ridx[i]] += r4;
12744324174eSBarry Smith       }
12754324174eSBarry Smith       b1 += bm4;
12764324174eSBarry Smith       b2 += bm4;
12774324174eSBarry Smith       b3 += bm4;
12784324174eSBarry Smith       b4 += bm4;
12794324174eSBarry Smith     }
12804324174eSBarry Smith     for (; col<cn; col++) {     /* over extra columns of C */
12814324174eSBarry Smith       colam = col*am;
12824324174eSBarry Smith       arm   = a->compressedrow.nrows;
12834324174eSBarry Smith       ii    = a->compressedrow.i;
12844324174eSBarry Smith       ridx  = a->compressedrow.rindex;
12854324174eSBarry Smith       for (i=0; i<arm; i++) {  /* over rows of C in those columns */
12864324174eSBarry Smith         r1 = 0.0;
12874324174eSBarry Smith         n  = ii[i+1] - ii[i];
12884324174eSBarry Smith         aj = a->j + ii[i];
12894324174eSBarry Smith         aa = a->a + ii[i];
12904324174eSBarry Smith 
12914324174eSBarry Smith         for (j=0; j<n; j++) {
12924324174eSBarry Smith           r1 += (*aa++)*b1[*aj++];
12934324174eSBarry Smith         }
1294a2ea699eSBarry Smith         c[colam + ridx[i]] += r1;
12954324174eSBarry Smith       }
12964324174eSBarry Smith       b1 += bm;
12974324174eSBarry Smith     }
12984324174eSBarry Smith   } else {
1299f32d5d43SBarry Smith     for (col=0; col<cn-4; col += 4) {  /* over columns of C */
1300f32d5d43SBarry Smith       colam = col*am;
1301f32d5d43SBarry Smith       for (i=0; i<am; i++) {        /* over rows of C in those columns */
1302f32d5d43SBarry Smith         r1 = r2 = r3 = r4 = 0.0;
1303f32d5d43SBarry Smith         n  = a->i[i+1] - a->i[i];
1304f32d5d43SBarry Smith         aj = a->j + a->i[i];
1305f32d5d43SBarry Smith         aa = a->a + a->i[i];
1306f32d5d43SBarry Smith         for (j=0; j<n; j++) {
1307f32d5d43SBarry Smith           r1 += (*aa)*b1[*aj];
1308f32d5d43SBarry Smith           r2 += (*aa)*b2[*aj];
1309f32d5d43SBarry Smith           r3 += (*aa)*b3[*aj];
1310f32d5d43SBarry Smith           r4 += (*aa++)*b4[*aj++];
1311f32d5d43SBarry Smith         }
1312f32d5d43SBarry Smith         c[colam + i]       += r1;
1313f32d5d43SBarry Smith         c[colam + am + i]  += r2;
1314f32d5d43SBarry Smith         c[colam + am2 + i] += r3;
1315f32d5d43SBarry Smith         c[colam + am3 + i] += r4;
1316f32d5d43SBarry Smith       }
1317f32d5d43SBarry Smith       b1 += bm4;
1318f32d5d43SBarry Smith       b2 += bm4;
1319f32d5d43SBarry Smith       b3 += bm4;
1320f32d5d43SBarry Smith       b4 += bm4;
1321f32d5d43SBarry Smith     }
1322f32d5d43SBarry Smith     for (; col<cn; col++) {     /* over extra columns of C */
1323a2ea699eSBarry Smith       colam = col*am;
1324f32d5d43SBarry Smith       for (i=0; i<am; i++) {  /* over rows of C in those columns */
1325f32d5d43SBarry Smith         r1 = 0.0;
1326f32d5d43SBarry Smith         n  = a->i[i+1] - a->i[i];
1327f32d5d43SBarry Smith         aj = a->j + a->i[i];
1328f32d5d43SBarry Smith         aa = a->a + a->i[i];
1329f32d5d43SBarry Smith 
1330f32d5d43SBarry Smith         for (j=0; j<n; j++) {
1331f32d5d43SBarry Smith           r1 += (*aa++)*b1[*aj++];
1332f32d5d43SBarry Smith         }
1333a2ea699eSBarry Smith         c[colam + i] += r1;
1334f32d5d43SBarry Smith       }
1335f32d5d43SBarry Smith       b1 += bm;
1336f32d5d43SBarry Smith     }
13374324174eSBarry Smith   }
1338dc0b31edSSatish Balay   ierr = PetscLogFlops(cn*2.0*a->nz);CHKERRQ(ierr);
13398c778c55SBarry Smith   ierr = MatDenseRestoreArray(C,&c);CHKERRQ(ierr);
13409123193aSHong Zhang   PetscFunctionReturn(0);
13419123193aSHong Zhang }
1342b1683b59SHong Zhang 
1343b9af6bddSHong Zhang PetscErrorCode  MatTransColoringApplySpToDen_SeqAIJ(MatTransposeColoring coloring,Mat B,Mat Btdense)
1344c8db22f6SHong Zhang {
1345c8db22f6SHong Zhang   PetscErrorCode ierr;
13462f5f1f90SHong Zhang   Mat_SeqAIJ     *b       = (Mat_SeqAIJ*)B->data;
13472f5f1f90SHong Zhang   Mat_SeqDense   *btdense = (Mat_SeqDense*)Btdense->data;
13482f5f1f90SHong Zhang   PetscInt       *bi      = b->i,*bj=b->j;
13492f5f1f90SHong Zhang   PetscInt       m        = Btdense->rmap->n,n=Btdense->cmap->n,j,k,l,col,anz,*btcol,brow,ncolumns;
13502f5f1f90SHong Zhang   MatScalar      *btval,*btval_den,*ba=b->a;
13512f5f1f90SHong Zhang   PetscInt       *columns=coloring->columns,*colorforcol=coloring->colorforcol,ncolors=coloring->ncolors;
1352c8db22f6SHong Zhang 
1353c8db22f6SHong Zhang   PetscFunctionBegin;
13542f5f1f90SHong Zhang   btval_den=btdense->v;
13552f5f1f90SHong Zhang   ierr     = PetscMemzero(btval_den,(m*n)*sizeof(MatScalar));CHKERRQ(ierr);
13562f5f1f90SHong Zhang   for (k=0; k<ncolors; k++) {
13572f5f1f90SHong Zhang     ncolumns = coloring->ncolumns[k];
13582f5f1f90SHong Zhang     for (l=0; l<ncolumns; l++) { /* insert a row of B to a column of Btdense */
1359d2853540SHong Zhang       col   = *(columns + colorforcol[k] + l);
13602f5f1f90SHong Zhang       btcol = bj + bi[col];
13612f5f1f90SHong Zhang       btval = ba + bi[col];
13622f5f1f90SHong Zhang       anz   = bi[col+1] - bi[col];
1363d2853540SHong Zhang       for (j=0; j<anz; j++) {
13642f5f1f90SHong Zhang         brow            = btcol[j];
13652f5f1f90SHong Zhang         btval_den[brow] = btval[j];
1366c8db22f6SHong Zhang       }
1367c8db22f6SHong Zhang     }
13682f5f1f90SHong Zhang     btval_den += m;
1369c8db22f6SHong Zhang   }
1370c8db22f6SHong Zhang   PetscFunctionReturn(0);
1371c8db22f6SHong Zhang }
1372c8db22f6SHong Zhang 
1373b9af6bddSHong Zhang PetscErrorCode MatTransColoringApplyDenToSp_SeqAIJ(MatTransposeColoring matcoloring,Mat Cden,Mat Csp)
13748972f759SHong Zhang {
13758972f759SHong Zhang   PetscErrorCode ierr;
1376b2d2b04fSHong Zhang   Mat_SeqAIJ     *csp = (Mat_SeqAIJ*)Csp->data;
1377077f23c2SHong Zhang   PetscScalar    *ca_den,*ca_den_ptr,*ca=csp->a;
1378f99a636bSHong Zhang   PetscInt       k,l,m=Cden->rmap->n,ncolors=matcoloring->ncolors;
1379e88777f2SHong Zhang   PetscInt       brows=matcoloring->brows,*den2sp=matcoloring->den2sp;
1380077f23c2SHong Zhang   PetscInt       nrows,*row,*idx;
1381077f23c2SHong Zhang   PetscInt       *rows=matcoloring->rows,*colorforrow=matcoloring->colorforrow;
13828972f759SHong Zhang 
13838972f759SHong Zhang   PetscFunctionBegin;
1384a3fe58edSHong Zhang   ierr   = MatDenseGetArray(Cden,&ca_den);CHKERRQ(ierr);
1385f99a636bSHong Zhang 
1386077f23c2SHong Zhang   if (brows > 0) {
1387077f23c2SHong Zhang     PetscInt *lstart,row_end,row_start;
1388980ae229SHong Zhang     lstart = matcoloring->lstart;
1389eeb4fd02SHong Zhang     ierr = PetscMemzero(lstart,ncolors*sizeof(PetscInt));CHKERRQ(ierr);
1390eeb4fd02SHong Zhang 
1391077f23c2SHong Zhang     row_end = brows;
1392eeb4fd02SHong Zhang     if (row_end > m) row_end = m;
1393077f23c2SHong Zhang     for (row_start=0; row_start<m; row_start+=brows) { /* loop over row blocks of Csp */
1394077f23c2SHong Zhang       ca_den_ptr = ca_den;
1395980ae229SHong Zhang       for (k=0; k<ncolors; k++) { /* loop over colors (columns of Cden) */
1396eeb4fd02SHong Zhang         nrows = matcoloring->nrows[k];
1397eeb4fd02SHong Zhang         row   = rows  + colorforrow[k];
1398eeb4fd02SHong Zhang         idx   = den2sp + colorforrow[k];
1399eeb4fd02SHong Zhang         for (l=lstart[k]; l<nrows; l++) {
1400eeb4fd02SHong Zhang           if (row[l] >= row_end) {
1401eeb4fd02SHong Zhang             lstart[k] = l;
1402eeb4fd02SHong Zhang             break;
1403eeb4fd02SHong Zhang           } else {
1404077f23c2SHong Zhang             ca[idx[l]] = ca_den_ptr[row[l]];
1405eeb4fd02SHong Zhang           }
1406eeb4fd02SHong Zhang         }
1407077f23c2SHong Zhang         ca_den_ptr += m;
1408eeb4fd02SHong Zhang       }
1409077f23c2SHong Zhang       row_end += brows;
1410eeb4fd02SHong Zhang       if (row_end > m) row_end = m;
1411eeb4fd02SHong Zhang     }
1412077f23c2SHong Zhang   } else { /* non-blocked impl: loop over columns of Csp - slow if Csp is large */
1413077f23c2SHong Zhang     ca_den_ptr = ca_den;
1414077f23c2SHong Zhang     for (k=0; k<ncolors; k++) {
1415077f23c2SHong Zhang       nrows = matcoloring->nrows[k];
1416077f23c2SHong Zhang       row   = rows  + colorforrow[k];
1417077f23c2SHong Zhang       idx   = den2sp + colorforrow[k];
1418077f23c2SHong Zhang       for (l=0; l<nrows; l++) {
1419077f23c2SHong Zhang         ca[idx[l]] = ca_den_ptr[row[l]];
1420077f23c2SHong Zhang       }
1421077f23c2SHong Zhang       ca_den_ptr += m;
1422077f23c2SHong Zhang     }
1423f99a636bSHong Zhang   }
1424f99a636bSHong Zhang 
1425a3fe58edSHong Zhang   ierr = MatDenseRestoreArray(Cden,&ca_den);CHKERRQ(ierr);
1426f99a636bSHong Zhang #if defined(PETSC_USE_INFO)
1427077f23c2SHong Zhang   if (matcoloring->brows > 0) {
1428f99a636bSHong Zhang     ierr = PetscInfo1(Csp,"Loop over %D row blocks for den2sp\n",brows);CHKERRQ(ierr);
1429e88777f2SHong Zhang   } else {
1430077f23c2SHong Zhang     ierr = PetscInfo(Csp,"Loop over colors/columns of Cden, inefficient for large sparse matrix product \n");CHKERRQ(ierr);
1431e88777f2SHong Zhang   }
1432f99a636bSHong Zhang #endif
14338972f759SHong Zhang   PetscFunctionReturn(0);
14348972f759SHong Zhang }
14358972f759SHong Zhang 
1436b9af6bddSHong Zhang PetscErrorCode MatTransposeColoringCreate_SeqAIJ(Mat mat,ISColoring iscoloring,MatTransposeColoring c)
1437b1683b59SHong Zhang {
1438b1683b59SHong Zhang   PetscErrorCode ierr;
1439e88777f2SHong Zhang   PetscInt       i,n,nrows,Nbs,j,k,m,ncols,col,cm;
14401a83f524SJed Brown   const PetscInt *is,*ci,*cj,*row_idx;
1441b28d80bdSHong Zhang   PetscInt       nis = iscoloring->n,*rowhit,bs = 1;
1442b1683b59SHong Zhang   IS             *isa;
1443b28d80bdSHong Zhang   Mat_SeqAIJ     *csp = (Mat_SeqAIJ*)mat->data;
1444e88777f2SHong Zhang   PetscInt       *colorforrow,*rows,*rows_i,*idxhit,*spidx,*den2sp,*den2sp_i;
1445e88777f2SHong Zhang   PetscInt       *colorforcol,*columns,*columns_i,brows;
1446e88777f2SHong Zhang   PetscBool      flg;
1447b1683b59SHong Zhang 
1448b1683b59SHong Zhang   PetscFunctionBegin;
1449b1683b59SHong Zhang   ierr = ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);CHKERRQ(ierr);
1450e99be685SHong Zhang 
14517ecccc15SHong Zhang   /* bs >1 is not being tested yet! */
1452e88777f2SHong Zhang   Nbs       = mat->cmap->N/bs;
1453b1683b59SHong Zhang   c->M      = mat->rmap->N/bs;  /* set total rows, columns and local rows */
1454e88777f2SHong Zhang   c->N      = Nbs;
1455e88777f2SHong Zhang   c->m      = c->M;
1456b1683b59SHong Zhang   c->rstart = 0;
1457e88777f2SHong Zhang   c->brows  = 100;
1458b1683b59SHong Zhang 
1459b1683b59SHong Zhang   c->ncolors = nis;
1460dcca6d9dSJed Brown   ierr = PetscMalloc3(nis,&c->ncolumns,nis,&c->nrows,nis+1,&colorforrow);CHKERRQ(ierr);
1461854ce69bSBarry Smith   ierr = PetscMalloc1(csp->nz+1,&rows);CHKERRQ(ierr);
1462854ce69bSBarry Smith   ierr = PetscMalloc1(csp->nz+1,&den2sp);CHKERRQ(ierr);
1463e88777f2SHong Zhang 
1464e88777f2SHong Zhang   brows = c->brows;
1465c5929fdfSBarry Smith   ierr = PetscOptionsGetInt(NULL,NULL,"-matden2sp_brows",&brows,&flg);CHKERRQ(ierr);
1466e88777f2SHong Zhang   if (flg) c->brows = brows;
1467eeb4fd02SHong Zhang   if (brows > 0) {
1468854ce69bSBarry Smith     ierr = PetscMalloc1(nis+1,&c->lstart);CHKERRQ(ierr);
1469e88777f2SHong Zhang   }
14702205254eSKarl Rupp 
1471d2853540SHong Zhang   colorforrow[0] = 0;
1472d2853540SHong Zhang   rows_i         = rows;
1473f99a636bSHong Zhang   den2sp_i       = den2sp;
1474b1683b59SHong Zhang 
1475854ce69bSBarry Smith   ierr = PetscMalloc1(nis+1,&colorforcol);CHKERRQ(ierr);
1476854ce69bSBarry Smith   ierr = PetscMalloc1(Nbs+1,&columns);CHKERRQ(ierr);
14772205254eSKarl Rupp 
1478d2853540SHong Zhang   colorforcol[0] = 0;
1479d2853540SHong Zhang   columns_i      = columns;
1480d2853540SHong Zhang 
1481077f23c2SHong Zhang   /* get column-wise storage of mat */
1482077f23c2SHong Zhang   ierr = MatGetColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr);
1483b1683b59SHong Zhang 
1484b28d80bdSHong Zhang   cm   = c->m;
1485854ce69bSBarry Smith   ierr = PetscMalloc1(cm+1,&rowhit);CHKERRQ(ierr);
1486854ce69bSBarry Smith   ierr = PetscMalloc1(cm+1,&idxhit);CHKERRQ(ierr);
1487f99a636bSHong Zhang   for (i=0; i<nis; i++) { /* loop over color */
1488b1683b59SHong Zhang     ierr = ISGetLocalSize(isa[i],&n);CHKERRQ(ierr);
1489b1683b59SHong Zhang     ierr = ISGetIndices(isa[i],&is);CHKERRQ(ierr);
14902205254eSKarl Rupp 
1491b1683b59SHong Zhang     c->ncolumns[i] = n;
1492b1683b59SHong Zhang     if (n) {
1493d2853540SHong Zhang       ierr = PetscMemcpy(columns_i,is,n*sizeof(PetscInt));CHKERRQ(ierr);
1494b1683b59SHong Zhang     }
1495d2853540SHong Zhang     colorforcol[i+1] = colorforcol[i] + n;
1496d2853540SHong Zhang     columns_i       += n;
1497b1683b59SHong Zhang 
1498b1683b59SHong Zhang     /* fast, crude version requires O(N*N) work */
1499e8354b3eSHong Zhang     ierr = PetscMemzero(rowhit,cm*sizeof(PetscInt));CHKERRQ(ierr);
1500e99be685SHong Zhang 
1501b7caf3d6SHong Zhang     for (j=0; j<n; j++) { /* loop over columns*/
1502b1683b59SHong Zhang       col     = is[j];
1503d2853540SHong Zhang       row_idx = cj + ci[col];
1504b1683b59SHong Zhang       m       = ci[col+1] - ci[col];
1505b7caf3d6SHong Zhang       for (k=0; k<m; k++) { /* loop over columns marking them in rowhit */
1506e99be685SHong Zhang         idxhit[*row_idx]   = spidx[ci[col] + k];
1507d2853540SHong Zhang         rowhit[*row_idx++] = col + 1;
1508b1683b59SHong Zhang       }
1509b1683b59SHong Zhang     }
1510b1683b59SHong Zhang     /* count the number of hits */
1511b1683b59SHong Zhang     nrows = 0;
1512e8354b3eSHong Zhang     for (j=0; j<cm; j++) {
1513b1683b59SHong Zhang       if (rowhit[j]) nrows++;
1514b1683b59SHong Zhang     }
1515b1683b59SHong Zhang     c->nrows[i]      = nrows;
1516d2853540SHong Zhang     colorforrow[i+1] = colorforrow[i] + nrows;
1517d2853540SHong Zhang 
1518b1683b59SHong Zhang     nrows = 0;
1519b7caf3d6SHong Zhang     for (j=0; j<cm; j++) { /* loop over rows */
1520b1683b59SHong Zhang       if (rowhit[j]) {
1521d2853540SHong Zhang         rows_i[nrows]   = j;
152212b89a43SHong Zhang         den2sp_i[nrows] = idxhit[j];
1523b1683b59SHong Zhang         nrows++;
1524b1683b59SHong Zhang       }
1525b1683b59SHong Zhang     }
1526e88777f2SHong Zhang     den2sp_i += nrows;
1527077f23c2SHong Zhang 
1528b1683b59SHong Zhang     ierr    = ISRestoreIndices(isa[i],&is);CHKERRQ(ierr);
1529f99a636bSHong Zhang     rows_i += nrows;
1530b1683b59SHong Zhang   }
15310298fd71SBarry Smith   ierr = MatRestoreColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr);
1532b28d80bdSHong Zhang   ierr = PetscFree(rowhit);CHKERRQ(ierr);
1533b1683b59SHong Zhang   ierr = ISColoringRestoreIS(iscoloring,&isa);CHKERRQ(ierr);
1534d2853540SHong Zhang   if (csp->nz != colorforrow[nis]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"csp->nz %d != colorforrow[nis] %d",csp->nz,colorforrow[nis]);
1535b28d80bdSHong Zhang 
1536d2853540SHong Zhang   c->colorforrow = colorforrow;
1537d2853540SHong Zhang   c->rows        = rows;
1538f99a636bSHong Zhang   c->den2sp      = den2sp;
1539d2853540SHong Zhang   c->colorforcol = colorforcol;
1540d2853540SHong Zhang   c->columns     = columns;
15412205254eSKarl Rupp 
1542f94ccd6cSHong Zhang   ierr = PetscFree(idxhit);CHKERRQ(ierr);
1543b1683b59SHong Zhang   PetscFunctionReturn(0);
1544b1683b59SHong Zhang }
1545d013fc79Sandi selinger 
1546d013fc79Sandi selinger /* Needed for MatMatMult_SeqAIJ_SeqAIJ_Combined() */
1547d013fc79Sandi selinger /* Append value to an array if the value is not present yet. A bitarray */
1548d013fc79Sandi selinger /* was used to determine if there is already an entry at this position. */
1549d013fc79Sandi selinger void appendToArray(PetscInt val, PetscInt *array, PetscInt *cnzi)
1550d013fc79Sandi selinger {
1551d013fc79Sandi selinger   array[(*cnzi)++] = val;
1552d013fc79Sandi selinger }
1553d013fc79Sandi selinger 
155473b67375Sandi selinger /* This algorithm combines the symbolic and numeric phase of matrix-matrix multiplication. */
1555d013fc79Sandi selinger PetscErrorCode MatMatMult_SeqAIJ_SeqAIJ_Combined(Mat A,Mat B,PetscReal fill,Mat *C)
1556d013fc79Sandi selinger {
1557d013fc79Sandi selinger   PetscErrorCode     ierr;
1558d013fc79Sandi selinger   PetscLogDouble     flops=0.0;
1559d013fc79Sandi selinger   Mat_SeqAIJ         *a  = (Mat_SeqAIJ*)A->data, *b = (Mat_SeqAIJ*)B->data, *c;
1560d013fc79Sandi selinger   const PetscInt     *ai = a->i,*bi = b->i, *aj = a->j;
1561d013fc79Sandi selinger   PetscInt           *ci,*cj,*cj_i;
1562d013fc79Sandi selinger   PetscScalar        *ca, *ca_i;
1563d013fc79Sandi selinger   PetscInt           c_maxmem = 0, a_maxrownnz = 0, a_rownnz, a_col;
1564d013fc79Sandi selinger   PetscInt           am = A->rmap->N, bn = B->cmap->N, bm = B->rmap->N;
1565d013fc79Sandi selinger   PetscInt           i, k, ndouble = 0;
1566d013fc79Sandi selinger   PetscReal          afill;
1567d013fc79Sandi selinger   PetscScalar        *c_row_val_dense;
1568d013fc79Sandi selinger   PetscBool          *c_row_idx_flags;
1569d013fc79Sandi selinger   PetscInt           *aj_i = a->j;
1570d013fc79Sandi selinger   PetscScalar        *aa_i = a->a;
1571d013fc79Sandi selinger 
1572d013fc79Sandi selinger   PetscFunctionBegin;
1573d013fc79Sandi selinger   /* Step 1: Determine upper bounds on memory for C */
157473b67375Sandi selinger   for (i=0; i<am; i++) { /* iterate over all rows of A */
1575d013fc79Sandi 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 */
1576d013fc79Sandi selinger     const PetscInt *acol = aj + ai[i]; /* column indices of nonzero entries in this row */
1577d013fc79Sandi selinger     a_rownnz = 0;
1578d013fc79Sandi selinger     for (k=0;k<anzi;++k) a_rownnz += bi[acol[k]+1] - bi[acol[k]];
1579d013fc79Sandi selinger     a_maxrownnz = PetscMax(a_maxrownnz, a_rownnz);
1580d013fc79Sandi selinger     c_maxmem += a_rownnz;
1581d013fc79Sandi selinger   }
1582d013fc79Sandi selinger   ierr = PetscMalloc1(am+1, &ci);               CHKERRQ(ierr);
1583d013fc79Sandi selinger   ierr = PetscMalloc1(bn, &c_row_val_dense);    CHKERRQ(ierr);
1584d013fc79Sandi selinger   ierr = PetscMalloc1(bn, &c_row_idx_flags);    CHKERRQ(ierr);
1585d013fc79Sandi selinger   ierr = PetscMalloc1(c_maxmem,&cj);            CHKERRQ(ierr);
1586d013fc79Sandi selinger   ierr = PetscMalloc1(c_maxmem,&ca);            CHKERRQ(ierr);
1587d013fc79Sandi selinger   ca_i = ca;
1588d013fc79Sandi selinger   cj_i = cj;
1589d013fc79Sandi selinger   ci[0] = 0;
159073b67375Sandi selinger   ierr = PetscMemzero(c_row_val_dense, bn * sizeof(PetscScalar));CHKERRQ(ierr);
159173b67375Sandi selinger   ierr = PetscMemzero(c_row_idx_flags, bn * sizeof(PetscBool));CHKERRQ(ierr);
1592d013fc79Sandi selinger   for (i=0; i<am; i++) {
1593d013fc79Sandi selinger     /* Step 2: Initialize the dense row vector for C  */
1594d013fc79Sandi 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 */
1595d013fc79Sandi selinger     PetscInt cnzi           = 0;
1596d013fc79Sandi selinger     PetscInt *bj_i;
1597d013fc79Sandi selinger     PetscScalar *ba_i;
1598d013fc79Sandi selinger 
1599d013fc79Sandi selinger     /* Step 3: Do the numerical calculations */
1600d013fc79Sandi selinger     for (a_col=0; a_col<anzi; a_col++) {          /* iterate over all non zero values in a row of A */
1601d013fc79Sandi selinger       PetscInt a_col_index = aj_i[a_col];
1602d013fc79Sandi selinger       const PetscInt bnzi = bi[a_col_index+1] - bi[a_col_index];
1603d013fc79Sandi selinger       flops += 2*bnzi;
1604d013fc79Sandi selinger       bj_i = b->j + bi[a_col_index];   /* points to the current row in bj */
1605d013fc79Sandi selinger       ba_i = b->a + bi[a_col_index];   /* points to the current row in ba */
1606d013fc79Sandi selinger       for (k=0; k<bnzi; ++k) { /* iterate over all non zeros of this row in B */
1607d013fc79Sandi selinger         if (c_row_idx_flags[ bj_i[k] ] == PETSC_FALSE) {
1608d013fc79Sandi selinger           appendToArray(bj_i[k], cj_i, &cnzi);
1609d013fc79Sandi selinger           c_row_idx_flags[ bj_i[k] ] = PETSC_TRUE;
1610d013fc79Sandi selinger         }
1611d013fc79Sandi selinger         c_row_val_dense[ bj_i[k] ] += aa_i[a_col] * ba_i[k];
1612d013fc79Sandi selinger       }
1613d013fc79Sandi selinger     }
1614d013fc79Sandi selinger 
1615d013fc79Sandi selinger     /* Sort array */
16163353a62bSKarl Rupp     ierr = PetscSortInt(cnzi, cj_i);CHKERRQ(ierr);
1617d013fc79Sandi selinger     /* Step 4 */
1618d013fc79Sandi selinger     for (k=0; k < cnzi; k++) {
1619d013fc79Sandi selinger       ca_i[k] = c_row_val_dense[cj_i[k]];
1620d013fc79Sandi selinger       c_row_val_dense[cj_i[k]] = 0.;
1621d013fc79Sandi selinger       c_row_idx_flags[cj_i[k]] = PETSC_FALSE;
1622d013fc79Sandi selinger     }
1623d013fc79Sandi selinger     /* terminate current row */
1624d013fc79Sandi selinger     aa_i += anzi;
1625d013fc79Sandi selinger     aj_i += anzi;
1626d013fc79Sandi selinger     ca_i += cnzi;
1627d013fc79Sandi selinger     cj_i += cnzi;
1628d013fc79Sandi selinger     ci[i+1] = ci[i] + cnzi;
1629d013fc79Sandi selinger     flops += cnzi;
1630d013fc79Sandi selinger   }
1631d013fc79Sandi selinger 
1632d013fc79Sandi selinger   /* Step 5 */
1633d013fc79Sandi selinger   /* Create the new matrix */
1634d013fc79Sandi selinger   ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,NULL,C);CHKERRQ(ierr);
1635d013fc79Sandi selinger   ierr = MatSetBlockSizesFromMats(*C,A,B);CHKERRQ(ierr);
1636*02fe1965SBarry Smith   ierr = MatSetType(*C,((PetscObject)A)->type_name);CHKERRQ(ierr);
1637d013fc79Sandi selinger 
1638d013fc79Sandi selinger   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
1639d013fc79Sandi selinger   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
1640d013fc79Sandi selinger   c          = (Mat_SeqAIJ*)((*C)->data);
1641d013fc79Sandi selinger   c->a       = ca;
1642d013fc79Sandi selinger   c->free_a  = PETSC_TRUE;
1643d013fc79Sandi selinger   c->free_ij = PETSC_TRUE;
1644d013fc79Sandi selinger   c->nonew   = 0;
1645d013fc79Sandi selinger 
1646d013fc79Sandi selinger   /* set MatInfo */
1647d013fc79Sandi selinger   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
1648d013fc79Sandi selinger   if (afill < 1.0) afill = 1.0;
1649d013fc79Sandi selinger   c->maxnz                     = ci[am];
1650d013fc79Sandi selinger   c->nz                        = ci[am];
1651d013fc79Sandi selinger   (*C)->info.mallocs           = ndouble;
1652d013fc79Sandi selinger   (*C)->info.fill_ratio_given  = fill;
1653d013fc79Sandi selinger   (*C)->info.fill_ratio_needed = afill;
1654d013fc79Sandi selinger 
165573b67375Sandi selinger   ierr = PetscFree(c_row_val_dense);CHKERRQ(ierr);
165673b67375Sandi selinger   ierr = PetscFree(c_row_idx_flags);CHKERRQ(ierr);
1657d013fc79Sandi selinger 
1658d013fc79Sandi selinger   ierr = MatAssemblyBegin(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1659d013fc79Sandi selinger   ierr = MatAssemblyEnd(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1660d013fc79Sandi selinger   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
1661d013fc79Sandi selinger   PetscFunctionReturn(0);
1662d013fc79Sandi selinger }
1663