xref: /petsc/src/mat/impls/aij/seq/matmatmult.c (revision 0b7e3e3d83be3c1111df233b75dab26d03d9912c)
1be1d678aSKris Buschelman 
2d50806bdSBarry Smith /*
32499ec78SHong Zhang   Defines matrix-matrix product routines for pairs of SeqAIJ matrices
4d50806bdSBarry Smith           C = A * B
5d50806bdSBarry Smith */
6d50806bdSBarry Smith 
7c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/aij.h> /*I "petscmat.h" I*/
8c6db04a5SJed Brown #include <../src/mat/utils/freespace.h>
9c6db04a5SJed Brown #include <petscbt.h>
10c6db04a5SJed Brown #include <../src/mat/impls/dense/seq/dense.h> /*I "petscmat.h" I*/
116284ec50SHong Zhang 
12ec01deb9SMatthew Knepley EXTERN_C_BEGIN
136284ec50SHong Zhang #undef __FUNCT__
145c66b693SKris Buschelman #define __FUNCT__ "MatMatMult_SeqAIJ_SeqAIJ"
1538baddfdSBarry Smith PetscErrorCode MatMatMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
1638baddfdSBarry Smith {
17dfbe8321SBarry Smith   PetscErrorCode ierr;
185c66b693SKris Buschelman 
195c66b693SKris Buschelman   PetscFunctionBegin;
2026be0446SHong Zhang   if (scall == MAT_INITIAL_MATRIX){
21519eb980SHong Zhang     /* ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr); */
2226be0446SHong Zhang     ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr);
23519eb980SHong Zhang     /* ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);   */
2426be0446SHong Zhang   }
25519eb980SHong Zhang   /* ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr); */
2626be0446SHong Zhang   ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr);
27519eb980SHong Zhang   /* ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr); */
285c66b693SKris Buschelman   PetscFunctionReturn(0);
295c66b693SKris Buschelman }
30ec01deb9SMatthew Knepley EXTERN_C_END
311c24bd37SHong Zhang 
3226be0446SHong Zhang #undef __FUNCT__
33b561aa9dSHong Zhang #define __FUNCT__ "MatGetSymbolicMatMatMult_SeqAIJ_SeqAIJ"
34b561aa9dSHong Zhang PetscErrorCode MatGetSymbolicMatMatMult_SeqAIJ_SeqAIJ(PetscInt am,PetscInt *Ai,PetscInt *Aj,PetscInt bm,PetscInt bn,PetscInt *Bi,PetscInt *Bj,PetscReal fill,PetscInt *Ci[],PetscInt *Cj[],PetscInt *nspacedouble)
3558c24d83SHong Zhang {
36dfbe8321SBarry Smith   PetscErrorCode     ierr;
37a1a86e44SBarry Smith   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
38b561aa9dSHong Zhang   PetscInt           *ai=Ai,*aj=Aj,*bi=Bi,*bj=Bj,*bjj,*ci,*cj;
39b561aa9dSHong Zhang   PetscInt           i,j,anzi,brow,bnzj,cnzi,nlnk,*lnk,ndouble=0;
40be0fcf8dSHong Zhang   PetscBT            lnkbt;
4158c24d83SHong Zhang 
4258c24d83SHong Zhang   PetscFunctionBegin;
4358c24d83SHong Zhang   /* Allocate ci array, arrays for fill computation and */
4458c24d83SHong Zhang   /* free space for accumulating nonzero column info */
4538baddfdSBarry Smith   ierr = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
4658c24d83SHong Zhang   ci[0] = 0;
4758c24d83SHong Zhang 
48be0fcf8dSHong Zhang   /* create and initialize a linked list */
49be0fcf8dSHong Zhang   nlnk = bn+1;
50be0fcf8dSHong Zhang   ierr = PetscLLCreate(bn,bn,nlnk,lnk,lnkbt);CHKERRQ(ierr);
5158c24d83SHong Zhang 
52c5db241fSHong Zhang   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
53a1a86e44SBarry Smith   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr);
5458c24d83SHong Zhang   current_space = free_space;
5558c24d83SHong Zhang 
5658c24d83SHong Zhang   /* Determine symbolic info for each row of the product: */
5758c24d83SHong Zhang   for (i=0;i<am;i++) {
5858c24d83SHong Zhang     anzi = ai[i+1] - ai[i];
5958c24d83SHong Zhang     cnzi = 0;
602d09714cSHong Zhang     j    = anzi;
61b561aa9dSHong Zhang     aj   = Aj + ai[i];
622d09714cSHong Zhang     while (j){/* assume cols are almost in increasing order, starting from its end saves computation */
632d09714cSHong Zhang       j--;
64b561aa9dSHong Zhang       brow = aj[j];
6558c24d83SHong Zhang       bnzj = bi[brow+1] - bi[brow];
6658c24d83SHong Zhang       bjj  = bj + bi[brow];
671c239cc6SHong Zhang       /* add non-zero cols of B into the sorted linked list lnk */
68be0fcf8dSHong Zhang       ierr = PetscLLAdd(bnzj,bjj,bn,nlnk,lnk,lnkbt);CHKERRQ(ierr);
691c239cc6SHong Zhang       cnzi += nlnk;
7058c24d83SHong Zhang     }
7158c24d83SHong Zhang 
7258c24d83SHong Zhang     /* If free space is not available, make more free space */
7358c24d83SHong Zhang     /* Double the amount of total space in the list */
7458c24d83SHong Zhang     if (current_space->local_remaining<cnzi) {
754238b7adSHong Zhang       ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,&current_space);CHKERRQ(ierr);
76b561aa9dSHong Zhang       ndouble++;
7758c24d83SHong Zhang     }
7858c24d83SHong Zhang 
79c5db241fSHong Zhang     /* Copy data into free space, then initialize lnk */
80be0fcf8dSHong Zhang     ierr = PetscLLClean(bn,bn,cnzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
81c5db241fSHong Zhang     current_space->array           += cnzi;
8258c24d83SHong Zhang     current_space->local_used      += cnzi;
8358c24d83SHong Zhang     current_space->local_remaining -= cnzi;
8458c24d83SHong Zhang     ci[i+1] = ci[i] + cnzi;
8558c24d83SHong Zhang   }
8658c24d83SHong Zhang 
8758c24d83SHong Zhang   /* Column indices are in the list of free space */
8858c24d83SHong Zhang   /* Allocate space for cj, initialize cj, and */
8958c24d83SHong Zhang   /* destroy list of free space and other temporary array(s) */
9038baddfdSBarry Smith   ierr = PetscMalloc((ci[am]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr);
91a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
92be0fcf8dSHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
9358c24d83SHong Zhang 
94b561aa9dSHong Zhang   *Ci           = ci;
95b561aa9dSHong Zhang   *Cj           = cj;
96b561aa9dSHong Zhang   *nspacedouble = ndouble;
97b561aa9dSHong Zhang   PetscFunctionReturn(0);
98b561aa9dSHong Zhang }
99b561aa9dSHong Zhang 
1000b7e3e3dSHong Zhang #define DENSEAXPY
101b561aa9dSHong Zhang #undef __FUNCT__
102b561aa9dSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ"
103b561aa9dSHong Zhang PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
104b561aa9dSHong Zhang {
105b561aa9dSHong Zhang   PetscErrorCode ierr;
106b561aa9dSHong Zhang   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
107b561aa9dSHong Zhang   PetscInt       *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*ci,*cj;
108b561aa9dSHong Zhang   PetscInt       am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N,nspacedouble;
109b561aa9dSHong Zhang   MatScalar      *ca;
110b561aa9dSHong Zhang   PetscReal      afill;
111b561aa9dSHong Zhang 
112b561aa9dSHong Zhang   PetscFunctionBegin;
113b561aa9dSHong Zhang   /* Get ci and cj */
114b561aa9dSHong Zhang   ierr = MatGetSymbolicMatMatMult_SeqAIJ_SeqAIJ(am,ai,aj,bm,bn,bi,bj,fill,&ci,&cj,&nspacedouble);CHKERRQ(ierr);
115b561aa9dSHong Zhang 
11658c24d83SHong Zhang   /* Allocate space for ca */
11758c24d83SHong Zhang   ierr = PetscMalloc((ci[am]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
11858c24d83SHong Zhang   ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr);
11958c24d83SHong Zhang 
12026be0446SHong Zhang   /* put together the new symbolic matrix */
1217adad957SLisandro Dalcin   ierr = MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,am,bn,ci,cj,ca,C);CHKERRQ(ierr);
12258c24d83SHong Zhang 
12358c24d83SHong Zhang   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
12458c24d83SHong Zhang   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
12558c24d83SHong Zhang   c = (Mat_SeqAIJ *)((*C)->data);
126e6b907acSBarry Smith   c->free_a   = PETSC_TRUE;
127e6b907acSBarry Smith   c->free_ij  = PETSC_TRUE;
12858c24d83SHong Zhang   c->nonew    = 0;
12958c24d83SHong Zhang 
1300b7e3e3dSHong Zhang #if defined(DENSEAXPY)
1310b7e3e3dSHong Zhang   ierr = PetscMalloc(B->cmap->N*sizeof(PetscScalar),&c->matmult_abdense);CHKERRQ(ierr);
1320b7e3e3dSHong Zhang   ierr = PetscMemzero(c->matmult_abdense,B->cmap->N*sizeof(PetscScalar));CHKERRQ(ierr);
1330b7e3e3dSHong Zhang #endif
1340b7e3e3dSHong Zhang 
1357212da7cSHong Zhang   /* set MatInfo */
1367212da7cSHong Zhang   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
137f2b054eeSHong Zhang   if (afill < 1.0) afill = 1.0;
1387212da7cSHong Zhang   c->maxnz                     = ci[am];
1397212da7cSHong Zhang   c->nz                        = ci[am];
1407212da7cSHong Zhang   (*C)->info.mallocs           = nspacedouble;
1417212da7cSHong Zhang   (*C)->info.fill_ratio_given  = fill;
1427212da7cSHong Zhang   (*C)->info.fill_ratio_needed = afill;
1437212da7cSHong Zhang 
1447212da7cSHong Zhang #if defined(PETSC_USE_INFO)
1457212da7cSHong Zhang   if (ci[am]) {
146f2b054eeSHong Zhang     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr);
147f2b054eeSHong Zhang     ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr);
148f2b054eeSHong Zhang   } else {
149f2b054eeSHong Zhang     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
150be0fcf8dSHong Zhang   }
151f2b054eeSHong Zhang #endif
15258c24d83SHong Zhang   PetscFunctionReturn(0);
15358c24d83SHong Zhang }
154d50806bdSBarry Smith 
15526be0446SHong Zhang #undef __FUNCT__
15626be0446SHong Zhang #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqAIJ"
157dfbe8321SBarry Smith PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
158d50806bdSBarry Smith {
159dfbe8321SBarry Smith   PetscErrorCode ierr;
160d13dce4bSSatish Balay   PetscLogDouble flops=0.0;
161d50806bdSBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ *)A->data;
162d50806bdSBarry Smith   Mat_SeqAIJ     *b = (Mat_SeqAIJ *)B->data;
163d50806bdSBarry Smith   Mat_SeqAIJ     *c = (Mat_SeqAIJ *)C->data;
16438baddfdSBarry Smith   PetscInt       *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j;
165d0f46423SBarry Smith   PetscInt       am=A->rmap->N,cm=C->rmap->N;
166519eb980SHong Zhang   PetscInt       i,j,k,anzi,bnzi,cnzi,brow;
167519eb980SHong Zhang   PetscScalar    *aa=a->a,*ba=b->a,*baj,*ca=c->a;
168519eb980SHong Zhang #if defined(DENSEAXPY)
1690b7e3e3dSHong Zhang   PetscScalar    *ab_dense=c->matmult_abdense;
170519eb980SHong Zhang #else
171519eb980SHong Zhang   PetscInt       nextb;
172519eb980SHong Zhang #endif
173d50806bdSBarry Smith 
174d50806bdSBarry Smith   PetscFunctionBegin;
175c124e916SHong Zhang   /* clean old values in C */
176c124e916SHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
177d50806bdSBarry Smith   /* Traverse A row-wise. */
178d50806bdSBarry Smith   /* Build the ith row in C by summing over nonzero columns in A, */
179d50806bdSBarry Smith   /* the rows of B corresponding to nonzeros of A. */
180519eb980SHong Zhang #if defined(DENSEAXPY)
181d50806bdSBarry Smith   for (i=0;i<am;i++) {
182d50806bdSBarry Smith     anzi = ai[i+1] - ai[i];
183519eb980SHong Zhang     cnzi = ci[i+1] - ci[i];
184d50806bdSBarry Smith     for (j=0;j<anzi;j++) {
185519eb980SHong Zhang       brow = aj[j];
186d50806bdSBarry Smith       bnzi = bi[brow+1] - bi[brow];
187d50806bdSBarry Smith       bjj  = bj + bi[brow];
188d50806bdSBarry Smith       baj  = ba + bi[brow];
189519eb980SHong Zhang       /* perform dense axpy */
190519eb980SHong Zhang       for (k=0; k<bnzi; k++) {
191519eb980SHong Zhang         ab_dense[bjj[k]] += aa[j]*baj[k];
192519eb980SHong Zhang       }
193519eb980SHong Zhang       flops += 2*bnzi;
194519eb980SHong Zhang     }
195519eb980SHong Zhang     for (k=0; k<cnzi; k++) {
196519eb980SHong Zhang       ca[k]          += ab_dense[cj[k]];
197519eb980SHong Zhang       ab_dense[cj[k]] = 0.0; /* zero ab_dense */
198519eb980SHong Zhang     }
199519eb980SHong Zhang     flops += cnzi;
200519eb980SHong Zhang     aj += anzi; aa += anzi;
201519eb980SHong Zhang     cj += cnzi; ca += cnzi;
202519eb980SHong Zhang   }
203519eb980SHong Zhang #else
204519eb980SHong Zhang   for (i=0;i<am;i++) {
205519eb980SHong Zhang     anzi = ai[i+1] - ai[i];
206519eb980SHong Zhang     cnzi = ci[i+1] - ci[i];
207519eb980SHong Zhang     for (j=0;j<anzi;j++) {
208519eb980SHong Zhang       brow = aj[j];
209519eb980SHong Zhang       bnzi = bi[brow+1] - bi[brow];
210519eb980SHong Zhang       bjj  = bj + bi[brow];
211519eb980SHong Zhang       baj  = ba + bi[brow];
212519eb980SHong Zhang       /* perform sparse axpy */
213c124e916SHong Zhang       nextb = 0;
214c124e916SHong Zhang       for (k=0; nextb<bnzi; k++) {
215c124e916SHong Zhang         if (cj[k] == bjj[nextb]){ /* ccol == bcol */
216519eb980SHong Zhang           ca[k] += aa[j]*baj[nextb++];
217c124e916SHong Zhang         }
218d50806bdSBarry Smith       }
219d50806bdSBarry Smith       flops += 2*bnzi;
220d50806bdSBarry Smith     }
221519eb980SHong Zhang     aj += anzi; aa += anzi;
222519eb980SHong Zhang     cj += cnzi; ca += cnzi;
223d50806bdSBarry Smith   }
224519eb980SHong Zhang #endif
225716bacf3SKris Buschelman   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
226716bacf3SKris Buschelman   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
227d50806bdSBarry Smith   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
228d50806bdSBarry Smith   PetscFunctionReturn(0);
229d50806bdSBarry Smith }
230bc011b1eSHong Zhang 
231d2853540SHong Zhang /* This routine is not used. Should be removed! */
232bc011b1eSHong Zhang #undef __FUNCT__
2336fc122caSHong Zhang #define __FUNCT__ "MatMatTransposeMult_SeqAIJ_SeqAIJ"
2346fc122caSHong Zhang PetscErrorCode MatMatTransposeMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
2355df89d91SHong Zhang {
236bc011b1eSHong Zhang   PetscErrorCode ierr;
237bc011b1eSHong Zhang 
238bc011b1eSHong Zhang   PetscFunctionBegin;
239bc011b1eSHong Zhang   if (scall == MAT_INITIAL_MATRIX){
2406fc122caSHong Zhang     ierr = MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr);
241bc011b1eSHong Zhang   }
2426fc122caSHong Zhang   ierr = MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr);
243bc011b1eSHong Zhang   PetscFunctionReturn(0);
244bc011b1eSHong Zhang }
245bc011b1eSHong Zhang 
246bc011b1eSHong Zhang #undef __FUNCT__
247e286af10SHong Zhang #define __FUNCT__ "PetscContainerDestroy_Mat_MatMatTransMult"
248e286af10SHong Zhang PetscErrorCode PetscContainerDestroy_Mat_MatMatTransMult(void *ptr)
2492128a86cSHong Zhang {
2502128a86cSHong Zhang   PetscErrorCode      ierr;
251e286af10SHong Zhang   Mat_MatMatTransMult *multtrans=(Mat_MatMatTransMult*)ptr;
2522128a86cSHong Zhang 
2532128a86cSHong Zhang   PetscFunctionBegin;
254b9af6bddSHong Zhang   ierr = MatTransposeColoringDestroy(&multtrans->matcoloring);CHKERRQ(ierr);
2552128a86cSHong Zhang   ierr = MatDestroy(&multtrans->Bt_den);CHKERRQ(ierr);
2562128a86cSHong Zhang   ierr = MatDestroy(&multtrans->ABt_den);CHKERRQ(ierr);
2572128a86cSHong Zhang   ierr = PetscFree(multtrans);CHKERRQ(ierr);
2582128a86cSHong Zhang   PetscFunctionReturn(0);
2592128a86cSHong Zhang }
2602128a86cSHong Zhang 
2612128a86cSHong Zhang #undef __FUNCT__
2622128a86cSHong Zhang #define __FUNCT__ "MatDestroy_SeqAIJ_MatMatMultTrans"
2632128a86cSHong Zhang PetscErrorCode MatDestroy_SeqAIJ_MatMatMultTrans(Mat A)
2642128a86cSHong Zhang {
2652128a86cSHong Zhang   PetscErrorCode      ierr;
2662128a86cSHong Zhang   PetscContainer      container;
267e286af10SHong Zhang   Mat_MatMatTransMult *multtrans=PETSC_NULL;
2682128a86cSHong Zhang 
2692128a86cSHong Zhang   PetscFunctionBegin;
270e286af10SHong Zhang   ierr = PetscObjectQuery((PetscObject)A,"Mat_MatMatTransMult",(PetscObject *)&container);CHKERRQ(ierr);
2712128a86cSHong Zhang   if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit");
2722128a86cSHong Zhang   ierr = PetscContainerGetPointer(container,(void **)&multtrans);CHKERRQ(ierr);
2732128a86cSHong Zhang   A->ops->destroy   = multtrans->destroy;
2742128a86cSHong Zhang   if (A->ops->destroy) {
2752128a86cSHong Zhang     ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
2762128a86cSHong Zhang   }
277e286af10SHong Zhang   ierr = PetscObjectCompose((PetscObject)A,"Mat_MatMatTransMult",0);CHKERRQ(ierr);
2782128a86cSHong Zhang   PetscFunctionReturn(0);
2792128a86cSHong Zhang }
2802128a86cSHong Zhang 
2812128a86cSHong Zhang #undef __FUNCT__
2826fc122caSHong Zhang #define __FUNCT__ "MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ"
2836fc122caSHong Zhang PetscErrorCode MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
284bc011b1eSHong Zhang {
2855df89d91SHong Zhang   PetscErrorCode      ierr;
28681d82fe4SHong Zhang   Mat                 Bt;
28781d82fe4SHong Zhang   PetscInt            *bti,*btj;
288e286af10SHong Zhang   Mat_MatMatTransMult *multtrans;
2892128a86cSHong Zhang   PetscContainer      container;
290d2853540SHong Zhang   PetscLogDouble      t0,tf,etime2=0.0;
291d2853540SHong Zhang 
29281d82fe4SHong Zhang   PetscFunctionBegin;
293d2853540SHong Zhang   ierr = PetscGetTime(&t0);CHKERRQ(ierr);
29481d82fe4SHong Zhang    /* create symbolic Bt */
29581d82fe4SHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(B,&bti,&btj);CHKERRQ(ierr);
29681d82fe4SHong Zhang   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,B->cmap->n,B->rmap->n,bti,btj,PETSC_NULL,&Bt);CHKERRQ(ierr);
29781d82fe4SHong Zhang 
29881d82fe4SHong Zhang   /* get symbolic C=A*Bt */
29981d82fe4SHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,Bt,fill,C);CHKERRQ(ierr);
30081d82fe4SHong Zhang 
3012128a86cSHong Zhang   /* create a supporting struct for reuse intermidiate dense matrices with matcoloring */
302e286af10SHong Zhang   ierr = PetscNew(Mat_MatMatTransMult,&multtrans);CHKERRQ(ierr);
3032128a86cSHong Zhang 
3042128a86cSHong Zhang   /* attach the supporting struct to C */
3052128a86cSHong Zhang   ierr = PetscContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr);
3062128a86cSHong Zhang   ierr = PetscContainerSetPointer(container,multtrans);CHKERRQ(ierr);
307e286af10SHong Zhang   ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_Mat_MatMatTransMult);CHKERRQ(ierr);
308e286af10SHong Zhang   ierr = PetscObjectCompose((PetscObject)(*C),"Mat_MatMatTransMult",(PetscObject)container);CHKERRQ(ierr);
3092128a86cSHong Zhang   ierr = PetscContainerDestroy(&container);CHKERRQ(ierr);
3102128a86cSHong Zhang 
3112128a86cSHong Zhang   multtrans->usecoloring = PETSC_FALSE;
3122128a86cSHong Zhang   multtrans->destroy = (*C)->ops->destroy;
3132128a86cSHong Zhang   (*C)->ops->destroy = MatDestroy_SeqAIJ_MatMatMultTrans;
3142128a86cSHong Zhang 
315d2853540SHong Zhang   ierr = PetscGetTime(&tf);CHKERRQ(ierr);
316d2853540SHong Zhang   etime2 += tf - t0;
317d2853540SHong Zhang 
318f400d4cbSHong Zhang   ierr = PetscOptionsGetBool(PETSC_NULL,"-matmattransmult_color",&multtrans->usecoloring,PETSC_NULL);CHKERRQ(ierr);
3192128a86cSHong Zhang   if (multtrans->usecoloring){
320b9af6bddSHong Zhang     /* Create MatTransposeColoring from symbolic C=A*B^T */
321b9af6bddSHong Zhang     MatTransposeColoring matcoloring;
3222128a86cSHong Zhang     ISColoring           iscoloring;
3232128a86cSHong Zhang     Mat                  Bt_dense,C_dense;
324d2853540SHong Zhang     PetscLogDouble       etime0=0.0,etime01=0.0,etime1=0.0;
325e8354b3eSHong Zhang 
326e8354b3eSHong Zhang     ierr = PetscGetTime(&t0);CHKERRQ(ierr);
327502de53fSHong Zhang     ierr = MatGetColoring(*C,MATCOLORINGLF,&iscoloring);CHKERRQ(ierr);
328d2853540SHong Zhang     ierr = PetscGetTime(&tf);CHKERRQ(ierr);
329d2853540SHong Zhang     etime0 += tf - t0;
330d2853540SHong Zhang 
331d2853540SHong Zhang     ierr = PetscGetTime(&t0);CHKERRQ(ierr);
332b9af6bddSHong Zhang     ierr = MatTransposeColoringCreate(*C,iscoloring,&matcoloring);CHKERRQ(ierr);
3332128a86cSHong Zhang     multtrans->matcoloring = matcoloring;
3342128a86cSHong Zhang     ierr = ISColoringDestroy(&iscoloring);CHKERRQ(ierr);
335e8354b3eSHong Zhang     ierr = PetscGetTime(&tf);CHKERRQ(ierr);
336d2853540SHong Zhang     etime01 += tf - t0;
3372128a86cSHong Zhang 
338e8354b3eSHong Zhang     ierr = PetscGetTime(&t0);CHKERRQ(ierr);
3392128a86cSHong Zhang     /* Create Bt_dense and C_dense = A*Bt_dense */
3402128a86cSHong Zhang     ierr = MatCreate(PETSC_COMM_SELF,&Bt_dense);CHKERRQ(ierr);
3412128a86cSHong Zhang     ierr = MatSetSizes(Bt_dense,A->cmap->n,matcoloring->ncolors,A->cmap->n,matcoloring->ncolors);CHKERRQ(ierr);
3422128a86cSHong Zhang     ierr = MatSetType(Bt_dense,MATSEQDENSE);CHKERRQ(ierr);
3432128a86cSHong Zhang     ierr = MatSeqDenseSetPreallocation(Bt_dense,PETSC_NULL);CHKERRQ(ierr);
3442128a86cSHong Zhang     Bt_dense->assembled = PETSC_TRUE;
3452128a86cSHong Zhang     multtrans->Bt_den = Bt_dense;
3462128a86cSHong Zhang 
3472128a86cSHong Zhang     ierr = MatCreate(PETSC_COMM_SELF,&C_dense);CHKERRQ(ierr);
3482128a86cSHong Zhang     ierr = MatSetSizes(C_dense,A->rmap->n,matcoloring->ncolors,A->rmap->n,matcoloring->ncolors);CHKERRQ(ierr);
3492128a86cSHong Zhang     ierr = MatSetType(C_dense,MATSEQDENSE);CHKERRQ(ierr);
3502128a86cSHong Zhang     ierr = MatSeqDenseSetPreallocation(C_dense,PETSC_NULL);CHKERRQ(ierr);
3512128a86cSHong Zhang     Bt_dense->assembled = PETSC_TRUE;
3522128a86cSHong Zhang     multtrans->ABt_den = C_dense;
353e8354b3eSHong Zhang     ierr = PetscGetTime(&tf);CHKERRQ(ierr);
354e8354b3eSHong Zhang     etime1 += tf - t0;
355f94ccd6cSHong Zhang 
356f94ccd6cSHong Zhang #if defined(PETSC_USE_INFO)
357f94ccd6cSHong Zhang     Mat_SeqAIJ *c=(Mat_SeqAIJ*)(*C)->data;
358f94ccd6cSHong Zhang     ierr = PetscInfo5(*C,"Bt_dense: %D,%D; Cnz %D / (cm*ncolors %D) = %g\n",A->cmap->n,matcoloring->ncolors,c->nz,A->rmap->n*matcoloring->ncolors,(PetscReal)(c->nz)/(A->rmap->n*matcoloring->ncolors));
359f94ccd6cSHong Zhang     ierr = PetscInfo5(*C,"Sym = GetColor %g + ColorCreate %g + MatDenseCreate %g + non-colorSym %g = %g\n",etime0,etime01,etime1,etime2,etime0+etime01+etime1+etime2);
360f94ccd6cSHong Zhang #endif
3612128a86cSHong Zhang   }
362*e99be685SHong Zhang   /* clean up */
363*e99be685SHong Zhang   ierr = MatDestroy(&Bt);CHKERRQ(ierr);
364*e99be685SHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(B,&bti,&btj);CHKERRQ(ierr);
3652128a86cSHong Zhang 
366f94ccd6cSHong Zhang 
367f94ccd6cSHong Zhang 
36881d82fe4SHong Zhang #if defined(INEFFICIENT_ALGORITHM)
36981d82fe4SHong Zhang   /* The algorithm below computes am*bm sparse inner-product - inefficient! It will be deleted later. */
3701fa1209cSHong Zhang   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
3711fa1209cSHong Zhang   Mat_SeqAIJ         *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
3721fa1209cSHong Zhang   PetscInt           *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*ci,*cj,*acol,*bcol;
3731fa1209cSHong Zhang   PetscInt           am=A->rmap->N,bm=B->rmap->N;
3741fa1209cSHong Zhang   PetscInt           i,j,anzi,bnzj,cnzi,nlnk,*lnk,nspacedouble=0,ka,kb,index[1];
3751fa1209cSHong Zhang   MatScalar          *ca;
3761fa1209cSHong Zhang   PetscBT            lnkbt;
3771fa1209cSHong Zhang   PetscReal          afill;
3785df89d91SHong Zhang 
3791fa1209cSHong Zhang   /* Allocate row pointer array ci  */
3801fa1209cSHong Zhang   ierr = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
3811fa1209cSHong Zhang   ci[0] = 0;
3821fa1209cSHong Zhang 
3831fa1209cSHong Zhang   /* Create and initialize a linked list for C columns */
3841fa1209cSHong Zhang   nlnk = bm+1;
3851fa1209cSHong Zhang   ierr = PetscLLCreate(bm,bm,nlnk,lnk,lnkbt);CHKERRQ(ierr);
3861fa1209cSHong Zhang 
3871fa1209cSHong Zhang   /* Initial FreeSpace with size fill*(nnz(A)+nnz(B)) */
3881fa1209cSHong Zhang   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr);
3891fa1209cSHong Zhang   current_space = free_space;
3901fa1209cSHong Zhang 
3911fa1209cSHong Zhang   /* Determine symbolic info for each row of the product A*B^T: */
3921fa1209cSHong Zhang   for (i=0; i<am; i++) {
3931fa1209cSHong Zhang     anzi = ai[i+1] - ai[i];
3941fa1209cSHong Zhang     cnzi = 0;
3951fa1209cSHong Zhang     acol = aj + ai[i];
3961fa1209cSHong Zhang     for (j=0; j<bm; j++){
3971fa1209cSHong Zhang       bnzj = bi[j+1] - bi[j];
3981fa1209cSHong Zhang       bcol= bj + bi[j];
39981d82fe4SHong Zhang       /* sparse inner-product c(i,j)=A[i,:]*B[j,:]^T */
4001fa1209cSHong Zhang       ka = 0; kb = 0;
4011fa1209cSHong Zhang       while (ka < anzi && kb < bnzj){
40281d82fe4SHong Zhang         while (acol[ka] < bcol[kb] && ka < anzi) ka++;
40381d82fe4SHong Zhang         if (ka == anzi) break;
40481d82fe4SHong Zhang         while (acol[ka] > bcol[kb] && kb < bnzj) kb++;
40581d82fe4SHong Zhang         if (kb == bnzj) break;
4061fa1209cSHong Zhang         if (acol[ka] == bcol[kb]){ /* add nonzero c(i,j) to lnk */
4071fa1209cSHong Zhang           index[0] = j;
4081fa1209cSHong Zhang           ierr = PetscLLAdd(1,index,bm,nlnk,lnk,lnkbt);CHKERRQ(ierr);
4091fa1209cSHong Zhang           cnzi++;
4101fa1209cSHong Zhang           break;
4111fa1209cSHong Zhang         }
4121fa1209cSHong Zhang       }
4131fa1209cSHong Zhang     }
4141fa1209cSHong Zhang 
4151fa1209cSHong Zhang     /* If free space is not available, make more free space */
4161fa1209cSHong Zhang     /* Double the amount of total space in the list */
4171fa1209cSHong Zhang     if (current_space->local_remaining<cnzi) {
4181fa1209cSHong Zhang       ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,&current_space);CHKERRQ(ierr);
4191fa1209cSHong Zhang       nspacedouble++;
4201fa1209cSHong Zhang     }
4211fa1209cSHong Zhang 
4221fa1209cSHong Zhang     /* Copy data into free space, then initialize lnk */
4231fa1209cSHong Zhang     ierr = PetscLLClean(bm,bm,cnzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
4241fa1209cSHong Zhang     current_space->array           += cnzi;
4251fa1209cSHong Zhang     current_space->local_used      += cnzi;
4261fa1209cSHong Zhang     current_space->local_remaining -= cnzi;
4271fa1209cSHong Zhang 
4281fa1209cSHong Zhang     ci[i+1] = ci[i] + cnzi;
4291fa1209cSHong Zhang   }
4301fa1209cSHong Zhang 
4311fa1209cSHong Zhang 
4321fa1209cSHong Zhang   /* Column indices are in the list of free space.
4331fa1209cSHong Zhang      Allocate array cj, copy column indices to cj, and destroy list of free space */
4341fa1209cSHong Zhang   ierr = PetscMalloc((ci[am]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr);
4351fa1209cSHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
4361fa1209cSHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
4371fa1209cSHong Zhang 
4381fa1209cSHong Zhang   /* Allocate space for ca */
4391fa1209cSHong Zhang   ierr = PetscMalloc((ci[am]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
4401fa1209cSHong Zhang   ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr);
4411fa1209cSHong Zhang 
4421fa1209cSHong Zhang   /* put together the new symbolic matrix */
4431fa1209cSHong Zhang   ierr = MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,am,bm,ci,cj,ca,C);CHKERRQ(ierr);
4441fa1209cSHong Zhang 
4451fa1209cSHong Zhang   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
4461fa1209cSHong Zhang   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
4471fa1209cSHong Zhang   c = (Mat_SeqAIJ *)((*C)->data);
4481fa1209cSHong Zhang   c->free_a   = PETSC_TRUE;
4491fa1209cSHong Zhang   c->free_ij  = PETSC_TRUE;
4501fa1209cSHong Zhang   c->nonew    = 0;
4511fa1209cSHong Zhang 
4521fa1209cSHong Zhang   /* set MatInfo */
4531fa1209cSHong Zhang   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
4541fa1209cSHong Zhang   if (afill < 1.0) afill = 1.0;
4551fa1209cSHong Zhang   c->maxnz                     = ci[am];
4561fa1209cSHong Zhang   c->nz                        = ci[am];
4571fa1209cSHong Zhang   (*C)->info.mallocs           = nspacedouble;
4581fa1209cSHong Zhang   (*C)->info.fill_ratio_given  = fill;
4591fa1209cSHong Zhang   (*C)->info.fill_ratio_needed = afill;
4601fa1209cSHong Zhang 
4611fa1209cSHong Zhang #if defined(PETSC_USE_INFO)
4621fa1209cSHong Zhang   if (ci[am]) {
4631fa1209cSHong Zhang     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr);
4646fc122caSHong Zhang     ierr = PetscInfo1((*C),"Use MatMatTransposeMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr);
4651fa1209cSHong Zhang   } else {
4661fa1209cSHong Zhang     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
4671fa1209cSHong Zhang   }
4681fa1209cSHong Zhang #endif
46981d82fe4SHong Zhang #endif
4705df89d91SHong Zhang   PetscFunctionReturn(0);
4715df89d91SHong Zhang }
4725df89d91SHong Zhang 
4736973769bSHong Zhang /* #define USE_ARRAY - for sparse dot product. Slower than !USE_ARRAY */
4745df89d91SHong Zhang #undef __FUNCT__
4756fc122caSHong Zhang #define __FUNCT__ "MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ"
4766fc122caSHong Zhang PetscErrorCode MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
4775df89d91SHong Zhang {
4785df89d91SHong Zhang   PetscErrorCode ierr;
4795df89d91SHong Zhang   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data;
480e2cac8caSJed Brown   PetscInt       *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,anzi,bnzj,nexta,nextb,*acol,*bcol,brow;
48181d82fe4SHong Zhang   PetscInt       cm=C->rmap->n,*ci=c->i,*cj=c->j,i,j,cnzi,*ccol;
4825df89d91SHong Zhang   PetscLogDouble flops=0.0;
4836973769bSHong Zhang   MatScalar      *aa=a->a,*aval,*ba=b->a,*bval,*ca=c->a,*cval;
484e286af10SHong Zhang   Mat_MatMatTransMult *multtrans;
4852128a86cSHong Zhang   PetscContainer      container;
4866973769bSHong Zhang #if defined(USE_ARRAY)
4876973769bSHong Zhang   MatScalar      *spdot;
4886973769bSHong Zhang #endif
4895df89d91SHong Zhang 
4905df89d91SHong Zhang   PetscFunctionBegin;
491e286af10SHong Zhang   ierr = PetscObjectQuery((PetscObject)C,"Mat_MatMatTransMult",(PetscObject *)&container);CHKERRQ(ierr);
4922128a86cSHong Zhang   if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit");
4932128a86cSHong Zhang   ierr  = PetscContainerGetPointer(container,(void **)&multtrans);CHKERRQ(ierr);
4942128a86cSHong Zhang   if (multtrans->usecoloring){
495b9af6bddSHong Zhang     MatTransposeColoring  matcoloring = multtrans->matcoloring;
496c8db22f6SHong Zhang     Mat                   Bt_dense;
497c8db22f6SHong Zhang     PetscInt              m,n;
498b2d2b04fSHong Zhang     PetscLogDouble t0,tf,etime0=0.0,etime1=0.0,etime2=0.0;
499a0b95ffbSSatish Balay     Mat C_dense = multtrans->ABt_den;
500a0b95ffbSSatish Balay 
5012128a86cSHong Zhang     Bt_dense = multtrans->Bt_den;
502c8db22f6SHong Zhang     ierr = MatGetLocalSize(Bt_dense,&m,&n);CHKERRQ(ierr);
503c8db22f6SHong Zhang 
504b9af6bddSHong Zhang     /* Get Bt_dense by Apply MatTransposeColoring to B */
5052128a86cSHong Zhang     ierr = PetscGetTime(&t0);CHKERRQ(ierr);
506b9af6bddSHong Zhang     ierr = MatTransColoringApplySpToDen(matcoloring,B,Bt_dense);CHKERRQ(ierr);
5072128a86cSHong Zhang     ierr = PetscGetTime(&tf);CHKERRQ(ierr);
508b2d2b04fSHong Zhang     etime0 += tf - t0;
509c8db22f6SHong Zhang 
510c8db22f6SHong Zhang     /* C_dense = A*Bt_dense */
5112128a86cSHong Zhang     ierr = PetscGetTime(&t0);CHKERRQ(ierr);
5122128a86cSHong Zhang     ierr = MatMatMultNumeric_SeqAIJ_SeqDense(A,Bt_dense,C_dense);CHKERRQ(ierr);
5132128a86cSHong Zhang     ierr = PetscGetTime(&tf);CHKERRQ(ierr);
5142128a86cSHong Zhang     etime2 += tf - t0;
515c8db22f6SHong Zhang 
5162128a86cSHong Zhang     /* Recover C from C_dense */
5172128a86cSHong Zhang     ierr = PetscGetTime(&t0);CHKERRQ(ierr);
518b9af6bddSHong Zhang     ierr = MatTransColoringApplyDenToSp(matcoloring,C_dense,C);CHKERRQ(ierr);
5192128a86cSHong Zhang     ierr = PetscGetTime(&tf);CHKERRQ(ierr);
5202128a86cSHong Zhang     etime1 += tf - t0;
521f94ccd6cSHong Zhang #if defined(PETSC_USE_INFO)
522f94ccd6cSHong Zhang     ierr = PetscInfo4(C,"Num = ColoringApply: %g %g + Mult_sp_dense %g = %g\n",etime0,etime1,etime2,etime0+etime1+etime2);
523f94ccd6cSHong Zhang #endif
524c8db22f6SHong Zhang     PetscFunctionReturn(0);
525c8db22f6SHong Zhang   }
526c8db22f6SHong Zhang 
5276973769bSHong Zhang #if defined(USE_ARRAY)
5286973769bSHong Zhang   /* allocate an array for implementing sparse inner-product */
529e2cac8caSJed Brown   ierr = PetscMalloc((A->cmap->n+1)*sizeof(MatScalar),&spdot);CHKERRQ(ierr);
530e2cac8caSJed Brown   ierr = PetscMemzero(spdot,(A->cmap->n+1)*sizeof(MatScalar));CHKERRQ(ierr);
5316973769bSHong Zhang #endif
5326973769bSHong Zhang 
53381d82fe4SHong Zhang   /* clear old values in C */
53481d82fe4SHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
5355df89d91SHong Zhang 
5361fa1209cSHong Zhang   for (i=0; i<cm; i++) {
53781d82fe4SHong Zhang     anzi = ai[i+1] - ai[i];
53881d82fe4SHong Zhang     acol = aj + ai[i];
5396973769bSHong Zhang     aval = aa + ai[i];
5401fa1209cSHong Zhang     cnzi = ci[i+1] - ci[i];
5411fa1209cSHong Zhang     ccol = cj + ci[i];
5426973769bSHong Zhang     cval = ca + ci[i];
5431fa1209cSHong Zhang     for (j=0; j<cnzi; j++){
54481d82fe4SHong Zhang       brow = ccol[j];
54581d82fe4SHong Zhang       bnzj = bi[brow+1] - bi[brow];
54681d82fe4SHong Zhang       bcol = bj + bi[brow];
5476973769bSHong Zhang       bval = ba + bi[brow];
5486973769bSHong Zhang 
54981d82fe4SHong Zhang       /* perform sparse inner-product c(i,j)=A[i,:]*B[j,:]^T */
5506973769bSHong Zhang #if defined(USE_ARRAY)
5516973769bSHong Zhang       /* put ba to spdot array */
5526973769bSHong Zhang       for (nextb=0; nextb<bnzj; nextb++) spdot[bcol[nextb]] = bval[nextb];
5536973769bSHong Zhang       /* c(i,j)=A[i,:]*B[j,:]^T */
5546973769bSHong Zhang       for (nexta=0; nexta<anzi; nexta++){
5556973769bSHong Zhang         cval[j] += spdot[acol[nexta]]*aval[nexta];
5566973769bSHong Zhang       }
5576973769bSHong Zhang       /* zero spdot array */
5586973769bSHong Zhang       for (nextb=0; nextb<bnzj; nextb++) spdot[bcol[nextb]] = 0.0;
5596973769bSHong Zhang #else
56081d82fe4SHong Zhang       nexta = 0; nextb = 0;
56181d82fe4SHong Zhang       while (nexta<anzi && nextb<bnzj){
56281d82fe4SHong Zhang         while (acol[nexta] < bcol[nextb] && nexta < anzi) nexta++;
56381d82fe4SHong Zhang         if (nexta == anzi) break;
56481d82fe4SHong Zhang         while (acol[nexta] > bcol[nextb] && nextb < bnzj) nextb++;
56581d82fe4SHong Zhang         if (nextb == bnzj) break;
56681d82fe4SHong Zhang         if (acol[nexta] == bcol[nextb]){
5676973769bSHong Zhang           cval[j] += aval[nexta]*bval[nextb];
56881d82fe4SHong Zhang           nexta++; nextb++;
56981d82fe4SHong Zhang           flops += 2;
5701fa1209cSHong Zhang         }
5711fa1209cSHong Zhang       }
5726973769bSHong Zhang #endif
57381d82fe4SHong Zhang     }
57481d82fe4SHong Zhang   }
57581d82fe4SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
57681d82fe4SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
57781d82fe4SHong Zhang   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
5786973769bSHong Zhang #if defined(USE_ARRAY)
5796973769bSHong Zhang   ierr = PetscFree(spdot);
5806973769bSHong Zhang #endif
5815df89d91SHong Zhang   PetscFunctionReturn(0);
5825df89d91SHong Zhang }
5835df89d91SHong Zhang 
5845df89d91SHong Zhang #undef __FUNCT__
58575648e8dSHong Zhang #define __FUNCT__ "MatTransposeMatMult_SeqAIJ_SeqAIJ"
58675648e8dSHong Zhang PetscErrorCode MatTransposeMatMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) {
5875df89d91SHong Zhang   PetscErrorCode ierr;
5885df89d91SHong Zhang 
5895df89d91SHong Zhang   PetscFunctionBegin;
5905df89d91SHong Zhang   if (scall == MAT_INITIAL_MATRIX){
59175648e8dSHong Zhang     ierr = MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr);
5925df89d91SHong Zhang   }
59375648e8dSHong Zhang   ierr = MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr);
5945df89d91SHong Zhang   PetscFunctionReturn(0);
5955df89d91SHong Zhang }
5965df89d91SHong Zhang 
5975df89d91SHong Zhang #undef __FUNCT__
59875648e8dSHong Zhang #define __FUNCT__ "MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ"
59975648e8dSHong Zhang PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
6005df89d91SHong Zhang {
601bc011b1eSHong Zhang   PetscErrorCode ierr;
602bc011b1eSHong Zhang   Mat            At;
60338baddfdSBarry Smith   PetscInt       *ati,*atj;
604bc011b1eSHong Zhang 
605bc011b1eSHong Zhang   PetscFunctionBegin;
606bc011b1eSHong Zhang   /* create symbolic At */
607bc011b1eSHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr);
608d0f46423SBarry Smith   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,A->cmap->n,A->rmap->n,ati,atj,PETSC_NULL,&At);CHKERRQ(ierr);
609bc011b1eSHong Zhang 
610bc011b1eSHong Zhang   /* get symbolic C=At*B */
611bc011b1eSHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(At,B,fill,C);CHKERRQ(ierr);
612bc011b1eSHong Zhang 
613bc011b1eSHong Zhang   /* clean up */
6146bf464f9SBarry Smith   ierr = MatDestroy(&At);CHKERRQ(ierr);
615bc011b1eSHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr);
616bc011b1eSHong Zhang   PetscFunctionReturn(0);
617bc011b1eSHong Zhang }
618bc011b1eSHong Zhang 
619bc011b1eSHong Zhang #undef __FUNCT__
62075648e8dSHong Zhang #define __FUNCT__ "MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ"
62175648e8dSHong Zhang PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
622bc011b1eSHong Zhang {
623bc011b1eSHong Zhang   PetscErrorCode ierr;
6240fbc74f4SHong Zhang   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data;
625d0f46423SBarry Smith   PetscInt       am=A->rmap->n,anzi,*ai=a->i,*aj=a->j,*bi=b->i,*bj,bnzi,nextb;
626d0f46423SBarry Smith   PetscInt       cm=C->rmap->n,*ci=c->i,*cj=c->j,crow,*cjj,i,j,k;
627d13dce4bSSatish Balay   PetscLogDouble flops=0.0;
6280fbc74f4SHong Zhang   MatScalar      *aa=a->a,*ba,*ca=c->a,*caj;
629bc011b1eSHong Zhang 
630bc011b1eSHong Zhang   PetscFunctionBegin;
631bc011b1eSHong Zhang   /* clear old values in C */
632bc011b1eSHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
633bc011b1eSHong Zhang 
634bc011b1eSHong Zhang   /* compute A^T*B using outer product (A^T)[:,i]*B[i,:] */
635bc011b1eSHong Zhang   for (i=0;i<am;i++) {
636bc011b1eSHong Zhang     bj   = b->j + bi[i];
637bc011b1eSHong Zhang     ba   = b->a + bi[i];
638bc011b1eSHong Zhang     bnzi = bi[i+1] - bi[i];
639bc011b1eSHong Zhang     anzi = ai[i+1] - ai[i];
640bc011b1eSHong Zhang     for (j=0; j<anzi; j++) {
641bc011b1eSHong Zhang       nextb = 0;
6420fbc74f4SHong Zhang       crow  = *aj++;
643bc011b1eSHong Zhang       cjj   = cj + ci[crow];
644bc011b1eSHong Zhang       caj   = ca + ci[crow];
645bc011b1eSHong Zhang       /* perform sparse axpy operation.  Note cjj includes bj. */
646bc011b1eSHong Zhang       for (k=0; nextb<bnzi; k++) {
6470fbc74f4SHong Zhang         if (cjj[k] == *(bj+nextb)) { /* ccol == bcol */
6480fbc74f4SHong Zhang           caj[k] += (*aa)*(*(ba+nextb));
649bc011b1eSHong Zhang           nextb++;
650bc011b1eSHong Zhang         }
651bc011b1eSHong Zhang       }
652bc011b1eSHong Zhang       flops += 2*bnzi;
6530fbc74f4SHong Zhang       aa++;
654bc011b1eSHong Zhang     }
655bc011b1eSHong Zhang   }
656bc011b1eSHong Zhang 
657bc011b1eSHong Zhang   /* Assemble the final matrix and clean up */
658bc011b1eSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
659bc011b1eSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
660bc011b1eSHong Zhang   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
661bc011b1eSHong Zhang   PetscFunctionReturn(0);
662bc011b1eSHong Zhang }
6639123193aSHong Zhang 
664ec01deb9SMatthew Knepley EXTERN_C_BEGIN
6659123193aSHong Zhang #undef __FUNCT__
6669123193aSHong Zhang #define __FUNCT__ "MatMatMult_SeqAIJ_SeqDense"
6679123193aSHong Zhang PetscErrorCode MatMatMult_SeqAIJ_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
6689123193aSHong Zhang {
6699123193aSHong Zhang   PetscErrorCode ierr;
6709123193aSHong Zhang 
6719123193aSHong Zhang   PetscFunctionBegin;
6729123193aSHong Zhang   if (scall == MAT_INITIAL_MATRIX){
6739123193aSHong Zhang     ierr = MatMatMultSymbolic_SeqAIJ_SeqDense(A,B,fill,C);CHKERRQ(ierr);
6749123193aSHong Zhang   }
6759123193aSHong Zhang   ierr = MatMatMultNumeric_SeqAIJ_SeqDense(A,B,*C);CHKERRQ(ierr);
6769123193aSHong Zhang   PetscFunctionReturn(0);
6779123193aSHong Zhang }
678ec01deb9SMatthew Knepley EXTERN_C_END
679edd81eecSMatthew Knepley 
6809123193aSHong Zhang #undef __FUNCT__
6819123193aSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqDense"
6829123193aSHong Zhang PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C)
6839123193aSHong Zhang {
6849123193aSHong Zhang   PetscErrorCode ierr;
6859123193aSHong Zhang 
6869123193aSHong Zhang   PetscFunctionBegin;
6875a586d82SBarry Smith   ierr = MatMatMultSymbolic_SeqDense_SeqDense(A,B,0.0,C);CHKERRQ(ierr);
6889123193aSHong Zhang   PetscFunctionReturn(0);
6899123193aSHong Zhang }
6909123193aSHong Zhang 
6919123193aSHong Zhang #undef __FUNCT__
6929123193aSHong Zhang #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqDense"
6939123193aSHong Zhang PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqDense(Mat A,Mat B,Mat C)
6949123193aSHong Zhang {
695f32d5d43SBarry Smith   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
6969123193aSHong Zhang   PetscErrorCode ierr;
697dd6ea824SBarry Smith   PetscScalar    *b,*c,r1,r2,r3,r4,*b1,*b2,*b3,*b4;
698dd6ea824SBarry Smith   MatScalar      *aa;
699d0f46423SBarry Smith   PetscInt       cm=C->rmap->n, cn=B->cmap->n, bm=B->rmap->n, col, i,j,n,*aj, am = A->rmap->n;
700f32d5d43SBarry Smith   PetscInt       am2 = 2*am, am3 = 3*am,  bm4 = 4*bm,colam;
7019123193aSHong Zhang 
7029123193aSHong Zhang   PetscFunctionBegin;
703f32d5d43SBarry Smith   if (!cm || !cn) PetscFunctionReturn(0);
704e32f2f54SBarry Smith   if (bm != A->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Number columns in A %D not equal rows in B %D\n",A->cmap->n,bm);
705e32f2f54SBarry 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);
706e32f2f54SBarry 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);
707f32d5d43SBarry Smith   ierr = MatGetArray(B,&b);CHKERRQ(ierr);
708f32d5d43SBarry Smith   ierr = MatGetArray(C,&c);CHKERRQ(ierr);
709f32d5d43SBarry Smith   b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm;
710f32d5d43SBarry Smith   for (col=0; col<cn-4; col += 4){  /* over columns of C */
711f32d5d43SBarry Smith     colam = col*am;
712f32d5d43SBarry Smith     for (i=0; i<am; i++) {        /* over rows of C in those columns */
713f32d5d43SBarry Smith       r1 = r2 = r3 = r4 = 0.0;
714f32d5d43SBarry Smith       n   = a->i[i+1] - a->i[i];
715f32d5d43SBarry Smith       aj  = a->j + a->i[i];
716f32d5d43SBarry Smith       aa  = a->a + a->i[i];
717f32d5d43SBarry Smith       for (j=0; j<n; j++) {
718f32d5d43SBarry Smith         r1 += (*aa)*b1[*aj];
719f32d5d43SBarry Smith         r2 += (*aa)*b2[*aj];
720f32d5d43SBarry Smith         r3 += (*aa)*b3[*aj];
721f32d5d43SBarry Smith         r4 += (*aa++)*b4[*aj++];
7229123193aSHong Zhang       }
723f32d5d43SBarry Smith       c[colam + i]       = r1;
724f32d5d43SBarry Smith       c[colam + am + i]  = r2;
725f32d5d43SBarry Smith       c[colam + am2 + i] = r3;
726f32d5d43SBarry Smith       c[colam + am3 + i] = r4;
727f32d5d43SBarry Smith     }
728f32d5d43SBarry Smith     b1 += bm4;
729f32d5d43SBarry Smith     b2 += bm4;
730f32d5d43SBarry Smith     b3 += bm4;
731f32d5d43SBarry Smith     b4 += bm4;
732f32d5d43SBarry Smith   }
733f32d5d43SBarry Smith   for (;col<cn; col++){     /* over extra columns of C */
734f32d5d43SBarry Smith     for (i=0; i<am; i++) {  /* over rows of C in those columns */
735f32d5d43SBarry Smith       r1 = 0.0;
736f32d5d43SBarry Smith       n   = a->i[i+1] - a->i[i];
737f32d5d43SBarry Smith       aj  = a->j + a->i[i];
738f32d5d43SBarry Smith       aa  = a->a + a->i[i];
739f32d5d43SBarry Smith 
740f32d5d43SBarry Smith       for (j=0; j<n; j++) {
741f32d5d43SBarry Smith         r1 += (*aa++)*b1[*aj++];
742f32d5d43SBarry Smith       }
743f32d5d43SBarry Smith       c[col*am + i]     = r1;
744f32d5d43SBarry Smith     }
745f32d5d43SBarry Smith     b1 += bm;
746f32d5d43SBarry Smith   }
747dc0b31edSSatish Balay   ierr = PetscLogFlops(cn*(2.0*a->nz));CHKERRQ(ierr);
748f32d5d43SBarry Smith   ierr = MatRestoreArray(B,&b);CHKERRQ(ierr);
749f32d5d43SBarry Smith   ierr = MatRestoreArray(C,&c);CHKERRQ(ierr);
750f32d5d43SBarry Smith   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
751f32d5d43SBarry Smith   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
752f32d5d43SBarry Smith   PetscFunctionReturn(0);
753f32d5d43SBarry Smith }
754f32d5d43SBarry Smith 
755f32d5d43SBarry Smith /*
7564324174eSBarry Smith    Note very similar to MatMult_SeqAIJ(), should generate both codes from same base
757f32d5d43SBarry Smith */
758f32d5d43SBarry Smith #undef __FUNCT__
759f32d5d43SBarry Smith #define __FUNCT__ "MatMatMultNumericAdd_SeqAIJ_SeqDense"
760f32d5d43SBarry Smith PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat A,Mat B,Mat C)
761f32d5d43SBarry Smith {
762f32d5d43SBarry Smith   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
763f32d5d43SBarry Smith   PetscErrorCode ierr;
764dd6ea824SBarry Smith   PetscScalar    *b,*c,r1,r2,r3,r4,*b1,*b2,*b3,*b4;
765dd6ea824SBarry Smith   MatScalar      *aa;
766d0f46423SBarry Smith   PetscInt       cm=C->rmap->n, cn=B->cmap->n, bm=B->rmap->n, col, i,j,n,*aj, am = A->rmap->n,*ii,arm;
7674324174eSBarry Smith   PetscInt       am2 = 2*am, am3 = 3*am,  bm4 = 4*bm,colam,*ridx;
768f32d5d43SBarry Smith 
769f32d5d43SBarry Smith   PetscFunctionBegin;
770f32d5d43SBarry Smith   if (!cm || !cn) PetscFunctionReturn(0);
771f32d5d43SBarry Smith   ierr = MatGetArray(B,&b);CHKERRQ(ierr);
772f32d5d43SBarry Smith   ierr = MatGetArray(C,&c);CHKERRQ(ierr);
773f32d5d43SBarry Smith   b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm;
7744324174eSBarry Smith 
7754324174eSBarry Smith   if (a->compressedrow.use){ /* use compressed row format */
7764324174eSBarry Smith     for (col=0; col<cn-4; col += 4){  /* over columns of C */
7774324174eSBarry Smith       colam = col*am;
7784324174eSBarry Smith       arm   = a->compressedrow.nrows;
7794324174eSBarry Smith       ii    = a->compressedrow.i;
7804324174eSBarry Smith       ridx  = a->compressedrow.rindex;
7814324174eSBarry Smith       for (i=0; i<arm; i++) {        /* over rows of C in those columns */
7824324174eSBarry Smith 	r1 = r2 = r3 = r4 = 0.0;
7834324174eSBarry Smith 	n   = ii[i+1] - ii[i];
7844324174eSBarry Smith 	aj  = a->j + ii[i];
7854324174eSBarry Smith 	aa  = a->a + ii[i];
7864324174eSBarry Smith 	for (j=0; j<n; j++) {
7874324174eSBarry Smith 	  r1 += (*aa)*b1[*aj];
7884324174eSBarry Smith 	  r2 += (*aa)*b2[*aj];
7894324174eSBarry Smith 	  r3 += (*aa)*b3[*aj];
7904324174eSBarry Smith 	  r4 += (*aa++)*b4[*aj++];
7914324174eSBarry Smith 	}
7924324174eSBarry Smith 	c[colam       + ridx[i]] += r1;
7934324174eSBarry Smith 	c[colam + am  + ridx[i]] += r2;
7944324174eSBarry Smith 	c[colam + am2 + ridx[i]] += r3;
7954324174eSBarry Smith 	c[colam + am3 + ridx[i]] += r4;
7964324174eSBarry Smith       }
7974324174eSBarry Smith       b1 += bm4;
7984324174eSBarry Smith       b2 += bm4;
7994324174eSBarry Smith       b3 += bm4;
8004324174eSBarry Smith       b4 += bm4;
8014324174eSBarry Smith     }
8024324174eSBarry Smith     for (;col<cn; col++){     /* over extra columns of C */
8034324174eSBarry Smith       colam = col*am;
8044324174eSBarry Smith       arm   = a->compressedrow.nrows;
8054324174eSBarry Smith       ii    = a->compressedrow.i;
8064324174eSBarry Smith       ridx  = a->compressedrow.rindex;
8074324174eSBarry Smith       for (i=0; i<arm; i++) {  /* over rows of C in those columns */
8084324174eSBarry Smith 	r1 = 0.0;
8094324174eSBarry Smith 	n   = ii[i+1] - ii[i];
8104324174eSBarry Smith 	aj  = a->j + ii[i];
8114324174eSBarry Smith 	aa  = a->a + ii[i];
8124324174eSBarry Smith 
8134324174eSBarry Smith 	for (j=0; j<n; j++) {
8144324174eSBarry Smith 	  r1 += (*aa++)*b1[*aj++];
8154324174eSBarry Smith 	}
8164324174eSBarry Smith 	c[col*am + ridx[i]] += r1;
8174324174eSBarry Smith       }
8184324174eSBarry Smith       b1 += bm;
8194324174eSBarry Smith     }
8204324174eSBarry Smith   } else {
821f32d5d43SBarry Smith     for (col=0; col<cn-4; col += 4){  /* over columns of C */
822f32d5d43SBarry Smith       colam = col*am;
823f32d5d43SBarry Smith       for (i=0; i<am; i++) {        /* over rows of C in those columns */
824f32d5d43SBarry Smith 	r1 = r2 = r3 = r4 = 0.0;
825f32d5d43SBarry Smith 	n   = a->i[i+1] - a->i[i];
826f32d5d43SBarry Smith 	aj  = a->j + a->i[i];
827f32d5d43SBarry Smith 	aa  = a->a + a->i[i];
828f32d5d43SBarry Smith 	for (j=0; j<n; j++) {
829f32d5d43SBarry Smith 	  r1 += (*aa)*b1[*aj];
830f32d5d43SBarry Smith 	  r2 += (*aa)*b2[*aj];
831f32d5d43SBarry Smith 	  r3 += (*aa)*b3[*aj];
832f32d5d43SBarry Smith 	  r4 += (*aa++)*b4[*aj++];
833f32d5d43SBarry Smith 	}
834f32d5d43SBarry Smith 	c[colam + i]       += r1;
835f32d5d43SBarry Smith 	c[colam + am + i]  += r2;
836f32d5d43SBarry Smith 	c[colam + am2 + i] += r3;
837f32d5d43SBarry Smith 	c[colam + am3 + i] += r4;
838f32d5d43SBarry Smith       }
839f32d5d43SBarry Smith       b1 += bm4;
840f32d5d43SBarry Smith       b2 += bm4;
841f32d5d43SBarry Smith       b3 += bm4;
842f32d5d43SBarry Smith       b4 += bm4;
843f32d5d43SBarry Smith     }
844f32d5d43SBarry Smith     for (;col<cn; col++){     /* over extra columns of C */
845f32d5d43SBarry Smith       for (i=0; i<am; i++) {  /* over rows of C in those columns */
846f32d5d43SBarry Smith 	r1 = 0.0;
847f32d5d43SBarry Smith 	n   = a->i[i+1] - a->i[i];
848f32d5d43SBarry Smith 	aj  = a->j + a->i[i];
849f32d5d43SBarry Smith 	aa  = a->a + a->i[i];
850f32d5d43SBarry Smith 
851f32d5d43SBarry Smith 	for (j=0; j<n; j++) {
852f32d5d43SBarry Smith 	  r1 += (*aa++)*b1[*aj++];
853f32d5d43SBarry Smith 	}
854f32d5d43SBarry Smith 	c[col*am + i]     += r1;
855f32d5d43SBarry Smith       }
856f32d5d43SBarry Smith       b1 += bm;
857f32d5d43SBarry Smith     }
8584324174eSBarry Smith   }
859dc0b31edSSatish Balay   ierr = PetscLogFlops(cn*2.0*a->nz);CHKERRQ(ierr);
860f32d5d43SBarry Smith   ierr = MatRestoreArray(B,&b);CHKERRQ(ierr);
861f32d5d43SBarry Smith   ierr = MatRestoreArray(C,&c);CHKERRQ(ierr);
8629123193aSHong Zhang   PetscFunctionReturn(0);
8639123193aSHong Zhang }
864b1683b59SHong Zhang 
865b1683b59SHong Zhang #undef __FUNCT__
866b9af6bddSHong Zhang #define __FUNCT__ "MatTransColoringApplySpToDen_SeqAIJ"
867b9af6bddSHong Zhang PetscErrorCode  MatTransColoringApplySpToDen_SeqAIJ(MatTransposeColoring coloring,Mat B,Mat Btdense)
868c8db22f6SHong Zhang {
869c8db22f6SHong Zhang   PetscErrorCode ierr;
8702f5f1f90SHong Zhang   Mat_SeqAIJ     *b = (Mat_SeqAIJ*)B->data;
8712f5f1f90SHong Zhang   Mat_SeqDense   *btdense = (Mat_SeqDense*)Btdense->data;
8722f5f1f90SHong Zhang   PetscInt       *bi=b->i,*bj=b->j;
8732f5f1f90SHong Zhang   PetscInt       m=Btdense->rmap->n,n=Btdense->cmap->n,j,k,l,col,anz,*btcol,brow,ncolumns;
8742f5f1f90SHong Zhang   MatScalar      *btval,*btval_den,*ba=b->a;
8752f5f1f90SHong Zhang   PetscInt       *columns=coloring->columns,*colorforcol=coloring->colorforcol,ncolors=coloring->ncolors;
876c8db22f6SHong Zhang 
877c8db22f6SHong Zhang   PetscFunctionBegin;
8782f5f1f90SHong Zhang   btval_den=btdense->v;
8792f5f1f90SHong Zhang   ierr = PetscMemzero(btval_den,(m*n)*sizeof(MatScalar));CHKERRQ(ierr);
8802f5f1f90SHong Zhang   for (k=0; k<ncolors; k++) {
8812f5f1f90SHong Zhang     ncolumns = coloring->ncolumns[k];
8822f5f1f90SHong Zhang     for (l=0; l<ncolumns; l++) { /* insert a row of B to a column of Btdense */
883d2853540SHong Zhang       col   = *(columns + colorforcol[k] + l);
8842f5f1f90SHong Zhang       btcol = bj + bi[col];
8852f5f1f90SHong Zhang       btval = ba + bi[col];
8862f5f1f90SHong Zhang       anz   = bi[col+1] - bi[col];
887d2853540SHong Zhang       for (j=0; j<anz; j++){
8882f5f1f90SHong Zhang         brow            = btcol[j];
8892f5f1f90SHong Zhang         btval_den[brow] = btval[j];
890c8db22f6SHong Zhang       }
891c8db22f6SHong Zhang     }
8922f5f1f90SHong Zhang     btval_den += m;
893c8db22f6SHong Zhang   }
894c8db22f6SHong Zhang   PetscFunctionReturn(0);
895c8db22f6SHong Zhang }
896c8db22f6SHong Zhang 
897c8db22f6SHong Zhang #undef __FUNCT__
898b9af6bddSHong Zhang #define __FUNCT__ "MatTransColoringApplyDenToSp_SeqAIJ"
899b9af6bddSHong Zhang PetscErrorCode MatTransColoringApplyDenToSp_SeqAIJ(MatTransposeColoring matcoloring,Mat Cden,Mat Csp)
9008972f759SHong Zhang {
9018972f759SHong Zhang   PetscErrorCode ierr;
902b2d2b04fSHong Zhang   Mat_SeqAIJ     *csp = (Mat_SeqAIJ*)Csp->data;
9032f5f1f90SHong Zhang   PetscInt       k,l,*row,*idx,m,ncolors=matcoloring->ncolors,nrows;
904b2d2b04fSHong Zhang   PetscScalar    *ca_den,*cp_den,*ca=csp->a;
9052f5f1f90SHong Zhang   PetscInt       *rows=matcoloring->rows,*spidx=matcoloring->columnsforspidx,*colorforrow=matcoloring->colorforrow;
9068972f759SHong Zhang 
9078972f759SHong Zhang   PetscFunctionBegin;
9088972f759SHong Zhang   ierr = MatGetLocalSize(Csp,&m,PETSC_NULL);CHKERRQ(ierr);
909b2d2b04fSHong Zhang   ierr = MatGetArray(Cden,&ca_den);CHKERRQ(ierr);
910b2d2b04fSHong Zhang   cp_den = ca_den;
9112f5f1f90SHong Zhang   for (k=0; k<ncolors; k++) {
9122f5f1f90SHong Zhang     nrows = matcoloring->nrows[k];
9132f5f1f90SHong Zhang     row   = rows  + colorforrow[k];
9142f5f1f90SHong Zhang     idx   = spidx + colorforrow[k];
9152f5f1f90SHong Zhang     for (l=0; l<nrows; l++){
9162f5f1f90SHong Zhang       ca[idx[l]] = cp_den[row[l]];
917b2d2b04fSHong Zhang     }
918b2d2b04fSHong Zhang     cp_den += m;
919b2d2b04fSHong Zhang   }
920b2d2b04fSHong Zhang   ierr = MatRestoreArray(Cden,&ca_den);CHKERRQ(ierr);
9218972f759SHong Zhang   PetscFunctionReturn(0);
9228972f759SHong Zhang }
9238972f759SHong Zhang 
924*e99be685SHong Zhang /*
925*e99be685SHong Zhang  MatGetColumnIJ_SeqAIJ_Color() and MatRestoreColumnIJ_SeqAIJ_Color() are customized from
926*e99be685SHong Zhang  MatGetColumnIJ_SeqAIJ() and MatRestoreColumnIJ_SeqAIJ() by adding an output
927*e99be685SHong Zhang  spidx[], index of a->j, to be used for setting 'columnsforspidx' in MatTransposeColoringCreate_SeqAIJ().
928*e99be685SHong Zhang  */
929*e99be685SHong Zhang #undef __FUNCT__
930*e99be685SHong Zhang #define __FUNCT__ "MatGetColumnIJ_SeqAIJ_Color"
931*e99be685SHong Zhang PetscErrorCode MatGetColumnIJ_SeqAIJ_Color(Mat A,PetscInt oshift,PetscBool  symmetric,PetscBool  inodecompressed,PetscInt *nn,PetscInt *ia[],PetscInt *ja[],PetscInt *spidx[],PetscBool  *done)
932*e99be685SHong Zhang {
933*e99be685SHong Zhang   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
934*e99be685SHong Zhang   PetscErrorCode ierr;
935*e99be685SHong Zhang   PetscInt       i,*collengths,*cia,*cja,n = A->cmap->n,m = A->rmap->n;
936*e99be685SHong Zhang   PetscInt       nz = a->i[m],row,*jj,mr,col;
937*e99be685SHong Zhang   PetscInt       *cspidx;
938*e99be685SHong Zhang 
939*e99be685SHong Zhang   PetscFunctionBegin;
940*e99be685SHong Zhang   *nn = n;
941*e99be685SHong Zhang   if (!ia) PetscFunctionReturn(0);
942*e99be685SHong Zhang   if (symmetric) {
943*e99be685SHong Zhang     SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"MatGetColumnIJ_SeqAIJ_Color() not supported for the case symmetric");
944*e99be685SHong Zhang     ierr = MatToSymmetricIJ_SeqAIJ(A->rmap->n,a->i,a->j,0,oshift,ia,ja);CHKERRQ(ierr);
945*e99be685SHong Zhang   } else {
946*e99be685SHong Zhang     ierr = PetscMalloc((n+1)*sizeof(PetscInt),&collengths);CHKERRQ(ierr);
947*e99be685SHong Zhang     ierr = PetscMemzero(collengths,n*sizeof(PetscInt));CHKERRQ(ierr);
948*e99be685SHong Zhang     ierr = PetscMalloc((n+1)*sizeof(PetscInt),&cia);CHKERRQ(ierr);
949*e99be685SHong Zhang     ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&cja);CHKERRQ(ierr);
950*e99be685SHong Zhang     ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&cspidx);CHKERRQ(ierr);
951*e99be685SHong Zhang     jj = a->j;
952*e99be685SHong Zhang     for (i=0; i<nz; i++) {
953*e99be685SHong Zhang       collengths[jj[i]]++;
954*e99be685SHong Zhang     }
955*e99be685SHong Zhang     cia[0] = oshift;
956*e99be685SHong Zhang     for (i=0; i<n; i++) {
957*e99be685SHong Zhang       cia[i+1] = cia[i] + collengths[i];
958*e99be685SHong Zhang     }
959*e99be685SHong Zhang     ierr = PetscMemzero(collengths,n*sizeof(PetscInt));CHKERRQ(ierr);
960*e99be685SHong Zhang     jj   = a->j;
961*e99be685SHong Zhang     for (row=0; row<m; row++) {
962*e99be685SHong Zhang       mr = a->i[row+1] - a->i[row];
963*e99be685SHong Zhang       for (i=0; i<mr; i++) {
964*e99be685SHong Zhang         col = *jj++;
965*e99be685SHong Zhang         cspidx[cia[col] + collengths[col] - oshift] = a->i[row] + i; /* index of a->j */
966*e99be685SHong Zhang         cja[cia[col] + collengths[col]++ - oshift] = row + oshift;
967*e99be685SHong Zhang       }
968*e99be685SHong Zhang     }
969*e99be685SHong Zhang     ierr = PetscFree(collengths);CHKERRQ(ierr);
970*e99be685SHong Zhang     *ia = cia; *ja = cja;
971*e99be685SHong Zhang     *spidx = cspidx;
972*e99be685SHong Zhang   }
973*e99be685SHong Zhang   PetscFunctionReturn(0);
974*e99be685SHong Zhang }
975*e99be685SHong Zhang 
976*e99be685SHong Zhang #undef __FUNCT__
977*e99be685SHong Zhang #define __FUNCT__ "MatRestoreColumnIJ_SeqAIJ_Color"
978*e99be685SHong Zhang PetscErrorCode MatRestoreColumnIJ_SeqAIJ_Color(Mat A,PetscInt oshift,PetscBool  symmetric,PetscBool  inodecompressed,PetscInt *n,PetscInt *ia[],PetscInt *ja[],PetscInt *spidx[],PetscBool  *done)
979*e99be685SHong Zhang {
980*e99be685SHong Zhang   PetscErrorCode ierr;
981*e99be685SHong Zhang 
982*e99be685SHong Zhang   PetscFunctionBegin;
983*e99be685SHong Zhang   if (!ia) PetscFunctionReturn(0);
984*e99be685SHong Zhang 
985*e99be685SHong Zhang   ierr = PetscFree(*ia);CHKERRQ(ierr);
986*e99be685SHong Zhang   ierr = PetscFree(*ja);CHKERRQ(ierr);
987*e99be685SHong Zhang   ierr = PetscFree(*spidx);CHKERRQ(ierr);
988*e99be685SHong Zhang   PetscFunctionReturn(0);
989*e99be685SHong Zhang }
990*e99be685SHong Zhang 
9918972f759SHong Zhang #undef __FUNCT__
992b9af6bddSHong Zhang #define __FUNCT__ "MatTransposeColoringCreate_SeqAIJ"
993b9af6bddSHong Zhang PetscErrorCode MatTransposeColoringCreate_SeqAIJ(Mat mat,ISColoring iscoloring,MatTransposeColoring c)
994b1683b59SHong Zhang {
995b1683b59SHong Zhang   PetscErrorCode ierr;
996d2853540SHong Zhang   PetscInt       i,n,nrows,N,j,k,m,*row_idx,*ci,*cj,ncols,col,cm;
997b1683b59SHong Zhang   const PetscInt *is;
998b28d80bdSHong Zhang   PetscInt       nis = iscoloring->n,*rowhit,bs = 1;
999b1683b59SHong Zhang   IS             *isa;
1000b1683b59SHong Zhang   PetscBool      done;
1001b1683b59SHong Zhang   PetscBool      flg1,flg2;
1002b28d80bdSHong Zhang   Mat_SeqAIJ     *csp = (Mat_SeqAIJ*)mat->data;
1003*e99be685SHong Zhang   PetscInt       *colorforrow,*rows,*rows_i,*columnsforspidx,*columnsforspidx_i,*idxhit,*spidx;
1004d2853540SHong Zhang   PetscInt       *colorforcol,*columns,*columns_i;
1005b1683b59SHong Zhang 
1006b1683b59SHong Zhang   PetscFunctionBegin;
1007b1683b59SHong Zhang   ierr = ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);CHKERRQ(ierr);
1008*e99be685SHong Zhang 
1009b1683b59SHong Zhang   /* this is ugly way to get blocksize but cannot call MatGetBlockSize() because AIJ can have bs > 1 */
1010b1683b59SHong Zhang   ierr = PetscTypeCompare((PetscObject)mat,MATSEQBAIJ,&flg1);CHKERRQ(ierr);
1011b1683b59SHong Zhang   ierr = PetscTypeCompare((PetscObject)mat,MATMPIBAIJ,&flg2);CHKERRQ(ierr);
1012b1683b59SHong Zhang   if (flg1 || flg2) {
1013b1683b59SHong Zhang     ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr);
1014b1683b59SHong Zhang   }
1015b1683b59SHong Zhang 
1016b1683b59SHong Zhang   N         = mat->cmap->N/bs;
1017b1683b59SHong Zhang   c->M      = mat->rmap->N/bs;  /* set total rows, columns and local rows */
1018b1683b59SHong Zhang   c->N      = mat->cmap->N/bs;
1019b1683b59SHong Zhang   c->m      = mat->rmap->N/bs;
1020b1683b59SHong Zhang   c->rstart = 0;
1021b1683b59SHong Zhang 
1022b1683b59SHong Zhang   c->ncolors = nis;
1023b1683b59SHong Zhang   ierr       = PetscMalloc(nis*sizeof(PetscInt),&c->ncolumns);CHKERRQ(ierr);
1024b28d80bdSHong Zhang   ierr       = PetscMalloc(nis*sizeof(PetscInt),&c->nrows);CHKERRQ(ierr);
1025d2853540SHong Zhang   ierr       = PetscMalloc2(csp->nz+1,PetscInt,&rows,csp->nz+1,PetscInt,&columnsforspidx);CHKERRQ(ierr);
1026d2853540SHong Zhang   ierr       = PetscMalloc((nis+1)*sizeof(PetscInt),&colorforrow);CHKERRQ(ierr);
1027d2853540SHong Zhang   colorforrow[0]    = 0;
1028d2853540SHong Zhang   rows_i            = rows;
1029d2853540SHong Zhang   columnsforspidx_i = columnsforspidx;
1030b1683b59SHong Zhang 
1031d2853540SHong Zhang   ierr = PetscMalloc((nis+1)*sizeof(PetscInt),&colorforcol);CHKERRQ(ierr);
1032d2853540SHong Zhang   ierr = PetscMalloc((N+1)*sizeof(PetscInt),&columns);CHKERRQ(ierr);
1033d2853540SHong Zhang   colorforcol[0] = 0;
1034d2853540SHong Zhang   columns_i      = columns;
1035d2853540SHong Zhang 
1036*e99be685SHong Zhang   ierr = MatGetColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,&done);CHKERRQ(ierr); /* column-wise storage of mat */
1037b1683b59SHong Zhang   if (!done) SETERRQ1(((PetscObject)mat)->comm,PETSC_ERR_SUP,"MatGetColumnIJ() not supported for matrix type %s",((PetscObject)mat)->type_name);
1038b1683b59SHong Zhang 
1039b28d80bdSHong Zhang   cm = c->m;
1040b28d80bdSHong Zhang   ierr = PetscMalloc((cm+1)*sizeof(PetscInt),&rowhit);CHKERRQ(ierr);
1041*e99be685SHong Zhang   ierr = PetscMalloc((cm+1)*sizeof(PetscInt),&idxhit);CHKERRQ(ierr);
1042b1683b59SHong Zhang   for (i=0; i<nis; i++) {
1043b1683b59SHong Zhang     ierr = ISGetLocalSize(isa[i],&n);CHKERRQ(ierr);
1044b1683b59SHong Zhang     ierr = ISGetIndices(isa[i],&is);CHKERRQ(ierr);
1045b1683b59SHong Zhang     c->ncolumns[i] = n;
1046b1683b59SHong Zhang     if (n) {
1047d2853540SHong Zhang       ierr = PetscMemcpy(columns_i,is,n*sizeof(PetscInt));CHKERRQ(ierr);
1048b1683b59SHong Zhang     }
1049d2853540SHong Zhang     colorforcol[i+1] = colorforcol[i] + n;
1050d2853540SHong Zhang     columns_i       += n;
1051b1683b59SHong Zhang 
1052b1683b59SHong Zhang     /* fast, crude version requires O(N*N) work */
1053e8354b3eSHong Zhang     ierr = PetscMemzero(rowhit,cm*sizeof(PetscInt));CHKERRQ(ierr);
1054*e99be685SHong Zhang 
1055b1683b59SHong Zhang     /* loop over columns*/
1056b1683b59SHong Zhang     for (j=0; j<n; j++) {
1057b1683b59SHong Zhang       col     = is[j];
1058d2853540SHong Zhang       row_idx = cj + ci[col];
1059b1683b59SHong Zhang       m       = ci[col+1] - ci[col];
1060b1683b59SHong Zhang       /* loop over columns marking them in rowhit */
1061b1683b59SHong Zhang       for (k=0; k<m; k++) {
1062*e99be685SHong Zhang         idxhit[*row_idx]   = spidx[ci[col] + k];
1063d2853540SHong Zhang         rowhit[*row_idx++] = col + 1;
1064b1683b59SHong Zhang       }
1065b1683b59SHong Zhang     }
1066b1683b59SHong Zhang     /* count the number of hits */
1067b1683b59SHong Zhang     nrows = 0;
1068e8354b3eSHong Zhang     for (j=0; j<cm; j++) {
1069b1683b59SHong Zhang       if (rowhit[j]) nrows++;
1070b1683b59SHong Zhang     }
1071b1683b59SHong Zhang     c->nrows[i]      = nrows;
1072d2853540SHong Zhang     colorforrow[i+1] = colorforrow[i] + nrows;
1073d2853540SHong Zhang 
1074b1683b59SHong Zhang     nrows       = 0;
1075e8354b3eSHong Zhang     for (j=0; j<cm; j++) {
1076b1683b59SHong Zhang       if (rowhit[j]) {
1077d2853540SHong Zhang         rows_i[nrows]            = j;
1078*e99be685SHong Zhang         columnsforspidx_i[nrows] = idxhit[j];
1079b1683b59SHong Zhang         nrows++;
1080b1683b59SHong Zhang       }
1081b1683b59SHong Zhang     }
1082b1683b59SHong Zhang     ierr = ISRestoreIndices(isa[i],&is);CHKERRQ(ierr);
1083d2853540SHong Zhang     rows_i += nrows; columnsforspidx_i += nrows;
1084b1683b59SHong Zhang   }
1085*e99be685SHong Zhang   ierr = MatRestoreColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,&done);CHKERRQ(ierr);
1086b28d80bdSHong Zhang   ierr = PetscFree(rowhit);CHKERRQ(ierr);
1087b1683b59SHong Zhang   ierr = ISColoringRestoreIS(iscoloring,&isa);CHKERRQ(ierr);
1088d2853540SHong Zhang #if defined(PETSC_USE_DEBUG)
1089d2853540SHong Zhang   if (csp->nz != colorforrow[nis]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"csp->nz %d != colorforrow[nis] %d",csp->nz,colorforrow[nis]);
1090d2853540SHong Zhang #endif
1091b28d80bdSHong Zhang 
1092d2853540SHong Zhang   c->colorforrow     = colorforrow;
1093d2853540SHong Zhang   c->rows            = rows;
1094d2853540SHong Zhang   c->columnsforspidx = columnsforspidx;
1095d2853540SHong Zhang   c->colorforcol     = colorforcol;
1096d2853540SHong Zhang   c->columns         = columns;
1097f94ccd6cSHong Zhang   ierr = PetscFree(idxhit);CHKERRQ(ierr);
1098b1683b59SHong Zhang   PetscFunctionReturn(0);
1099b1683b59SHong Zhang }
1100