xref: /petsc/src/mat/impls/aij/seq/matmatmult.c (revision f747e1ae54eac4c971caa6b7e82576b5c3be8258)
1d50806bdSBarry Smith /*
2a50f8bf6SHong Zhang   Defines matrix-matrix product routines for pairs of AIJ matrices
3d50806bdSBarry Smith           C = A * B
4d50806bdSBarry Smith */
5d50806bdSBarry Smith 
6c1f4806aSKris Buschelman #include "src/mat/impls/aij/seq/aij.h" /*I "petscmat.h" I*/
770f19b1fSKris Buschelman #include "src/mat/utils/freespace.h"
82d09714cSHong Zhang #include "src/mat/impls/aij/mpi/mpiaij.h"
9d50806bdSBarry Smith 
1026be0446SHong Zhang typedef struct { /* used by MatMatMult_MPIAIJ_MPIAIJ for reusing symbolic mat product */
1126be0446SHong Zhang   IS     isrowa,isrowb,iscolb;
1226be0446SHong Zhang   Mat    *aseq,*bseq,C_seq;
1326be0446SHong Zhang } Mat_MatMatMultMPI;
1426be0446SHong Zhang 
152216b3a4SKris Buschelman static int logkey_matmatmult_symbolic = 0;
162216b3a4SKris Buschelman static int logkey_matmatmult_numeric  = 0;
172216b3a4SKris Buschelman 
18c1f4806aSKris Buschelman #undef __FUNCT__
19c1f4806aSKris Buschelman #define __FUNCT__ "MatMatMult"
205c66b693SKris Buschelman /*@
215c66b693SKris Buschelman    MatMatMult - Performs Matrix-Matrix Multiplication C=A*B.
2294e3eecaSKris Buschelman 
235c66b693SKris Buschelman    Collective on Mat
24d50806bdSBarry Smith 
255c66b693SKris Buschelman    Input Parameters:
265c66b693SKris Buschelman +  A - the left matrix
271c741599SHong Zhang .  B - the right matrix
281c741599SHong Zhang .  scall - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX
29c5db241fSHong Zhang -  fill - expected fill as ratio of nnz(C)/(nnz(A) + nnz(B))
305c66b693SKris Buschelman 
315c66b693SKris Buschelman    Output Parameters:
325c66b693SKris Buschelman .  C - the product matrix
335c66b693SKris Buschelman 
345c66b693SKris Buschelman    Notes:
355c66b693SKris Buschelman    C will be created and must be destroyed by the user with MatDestroy().
365c66b693SKris Buschelman 
374d3841fdSKris Buschelman    This routine is currently only implemented for pairs of SeqAIJ matrices and classes
384d3841fdSKris Buschelman    which inherit from SeqAIJ.  C will be of type MATSEQAIJ.
395c66b693SKris Buschelman 
405c66b693SKris Buschelman    Level: intermediate
415c66b693SKris Buschelman 
425c66b693SKris Buschelman .seealso: MatMatMultSymbolic(),MatMatMultNumeric()
435c66b693SKris Buschelman @*/
441c741599SHong Zhang int MatMatMult(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
452d09714cSHong Zhang {
46d50806bdSBarry Smith   int  ierr;
472d09714cSHong Zhang 
48d50806bdSBarry Smith   PetscFunctionBegin;
494482741eSBarry Smith   PetscValidHeaderSpecific(A,MAT_COOKIE,1);
50c9780b6fSBarry Smith   PetscValidType(A,1);
515c66b693SKris Buschelman   MatPreallocated(A);
525c66b693SKris Buschelman   if (!A->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
535c66b693SKris Buschelman   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
544482741eSBarry Smith   PetscValidHeaderSpecific(B,MAT_COOKIE,2);
55c9780b6fSBarry Smith   PetscValidType(B,2);
565c66b693SKris Buschelman   MatPreallocated(B);
575c66b693SKris Buschelman   if (!B->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
585c66b693SKris Buschelman   if (B->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
594482741eSBarry Smith   PetscValidPointer(C,3);
605c66b693SKris Buschelman   if (B->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",B->M,A->N);
61d50806bdSBarry Smith 
621c24bd37SHong Zhang   if (fill <=0.0) SETERRQ1(PETSC_ERR_ARG_SIZ,"fill=%g must be > 0.0",fill);
63c5db241fSHong Zhang 
646284ec50SHong Zhang   ierr = PetscLogEventBegin(MAT_MatMult,A,B,0,0);CHKERRQ(ierr);
651c741599SHong Zhang   ierr = (*A->ops->matmult)(A,B,scall,fill,C);CHKERRQ(ierr);
666284ec50SHong Zhang   ierr = PetscLogEventEnd(MAT_MatMult,A,B,0,0);CHKERRQ(ierr);
674d3841fdSKris Buschelman 
686284ec50SHong Zhang   PetscFunctionReturn(0);
696284ec50SHong Zhang }
706284ec50SHong Zhang 
716284ec50SHong Zhang #undef __FUNCT__
726284ec50SHong Zhang #define __FUNCT__ "MatMatMult_MPIAIJ_MPIAIJ"
731c741599SHong Zhang int MatMatMult_MPIAIJ_MPIAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill, Mat *C)
742d09714cSHong Zhang {
7526be0446SHong Zhang   int ierr;
766284ec50SHong Zhang 
776284ec50SHong Zhang   PetscFunctionBegin;
7826be0446SHong Zhang   if (scall == MAT_INITIAL_MATRIX){
79d6bb3c2dSHong Zhang     ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ(A,B,fill,C);CHKERRQ(ierr);/* numeric product is computed as well */
80d6bb3c2dSHong Zhang   } else if (scall == MAT_REUSE_MATRIX){
8126be0446SHong Zhang     ierr = MatMatMultNumeric_MPIAIJ_MPIAIJ(A,B,*C);CHKERRQ(ierr);
82d6bb3c2dSHong Zhang   } else {
83d6bb3c2dSHong Zhang     SETERRQ1(PETSC_ERR_ARG_WRONG,"Invalid MatReuse %d",scall);
84d6bb3c2dSHong Zhang   }
85d50806bdSBarry Smith   PetscFunctionReturn(0);
86d50806bdSBarry Smith }
875c66b693SKris Buschelman 
885c66b693SKris Buschelman #undef __FUNCT__
895c66b693SKris Buschelman #define __FUNCT__ "MatMatMult_SeqAIJ_SeqAIJ"
901c741599SHong Zhang int MatMatMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) {
915c66b693SKris Buschelman   int ierr;
925c66b693SKris Buschelman 
935c66b693SKris Buschelman   PetscFunctionBegin;
9426be0446SHong Zhang   if (scall == MAT_INITIAL_MATRIX){
9526be0446SHong Zhang     ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr);
9626be0446SHong Zhang   }
9726be0446SHong Zhang   ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr);
985c66b693SKris Buschelman   PetscFunctionReturn(0);
995c66b693SKris Buschelman }
1005c66b693SKris Buschelman 
101c1f4806aSKris Buschelman #undef __FUNCT__
102c1f4806aSKris Buschelman #define __FUNCT__ "MatMatMultSymbolic"
1035c66b693SKris Buschelman /*@
1045c66b693SKris Buschelman    MatMatMultSymbolic - Performs construction, preallocation, and computes the ij structure
1055c66b693SKris Buschelman    of the matrix-matrix product C=A*B.  Call this routine before calling MatMatMultNumeric().
1065c66b693SKris Buschelman 
1075c66b693SKris Buschelman    Collective on Mat
1085c66b693SKris Buschelman 
1095c66b693SKris Buschelman    Input Parameters:
1105c66b693SKris Buschelman +  A - the left matrix
111c5db241fSHong Zhang .  B - the right matrix
112c5db241fSHong Zhang -  fill - expected fill as ratio of nnz(C)/(nnz(A) + nnz(B))
1135c66b693SKris Buschelman 
1145c66b693SKris Buschelman    Output Parameters:
1155c66b693SKris Buschelman .  C - the matrix containing the ij structure of product matrix
1165c66b693SKris Buschelman 
1175c66b693SKris Buschelman    Notes:
1184d3841fdSKris Buschelman    C will be created as a MATSEQAIJ matrix and must be destroyed by the user with MatDestroy().
1195c66b693SKris Buschelman 
1204d3841fdSKris Buschelman    This routine is currently only implemented for SeqAIJ matrices and classes which inherit from SeqAIJ.
1215c66b693SKris Buschelman 
1225c66b693SKris Buschelman    Level: intermediate
1235c66b693SKris Buschelman 
1245c66b693SKris Buschelman .seealso: MatMatMult(),MatMatMultNumeric()
1255c66b693SKris Buschelman @*/
126c5db241fSHong Zhang int MatMatMultSymbolic(Mat A,Mat B,PetscReal fill,Mat *C) {
1275c66b693SKris Buschelman   /* Perhaps this "interface" routine should be moved into the interface directory.*/
1285c66b693SKris Buschelman   /* To facilitate implementations with varying types, QueryFunction is used.*/
1295c66b693SKris Buschelman   /* It is assumed that implementations will be composed as "MatMatMultSymbolic_<type of A><type of B>". */
1305c66b693SKris Buschelman   int  ierr;
1314d3841fdSKris Buschelman   char symfunct[80];
132c5db241fSHong Zhang   int  (*symbolic)(Mat,Mat,PetscReal,Mat *);
1335c66b693SKris Buschelman 
1345c66b693SKris Buschelman   PetscFunctionBegin;
1354482741eSBarry Smith   PetscValidHeaderSpecific(A,MAT_COOKIE,1);
136c9780b6fSBarry Smith   PetscValidType(A,1);
1375c66b693SKris Buschelman   MatPreallocated(A);
1385c66b693SKris Buschelman   if (!A->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
1395c66b693SKris Buschelman   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1405c66b693SKris Buschelman 
1414482741eSBarry Smith   PetscValidHeaderSpecific(B,MAT_COOKIE,2);
142c9780b6fSBarry Smith   PetscValidType(B,2);
1435c66b693SKris Buschelman   MatPreallocated(B);
1445c66b693SKris Buschelman   if (!B->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
1455c66b693SKris Buschelman   if (B->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1464482741eSBarry Smith   PetscValidPointer(C,3);
1474482741eSBarry Smith 
1485c66b693SKris Buschelman   if (B->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",B->M,A->N);
1491c24bd37SHong Zhang   if (fill <=0.0) SETERRQ1(PETSC_ERR_ARG_SIZ,"fill=%g must be > 0.0",fill);
1505c66b693SKris Buschelman 
1514d3841fdSKris Buschelman   /* Currently only _seqaijseqaij is implemented, so just query for it in A and B. */
1524d3841fdSKris Buschelman   /* When other implementations exist, attack the multiple dispatch problem. */
1534d3841fdSKris Buschelman   ierr = PetscStrcpy(symfunct,"MatMatMultSymbolic_seqaijseqaij");CHKERRQ(ierr);
1544d3841fdSKris Buschelman   ierr = PetscObjectQueryFunction((PetscObject)B,symfunct,(PetscVoidFunction)&symbolic);CHKERRQ(ierr);
1554d3841fdSKris Buschelman   if (!symbolic) SETERRQ1(PETSC_ERR_SUP,"C=A*B not implemented for B of type %s",B->type_name);
1564d3841fdSKris Buschelman   ierr = PetscObjectQueryFunction((PetscObject)A,symfunct,(PetscVoidFunction)&symbolic);CHKERRQ(ierr);
1574d3841fdSKris Buschelman   if (!symbolic) SETERRQ1(PETSC_ERR_SUP,"C=A*B not implemented for A of type %s",A->type_name);
158c5db241fSHong Zhang   ierr = (*symbolic)(A,B,fill,C);CHKERRQ(ierr);
1595c66b693SKris Buschelman 
1605c66b693SKris Buschelman   PetscFunctionReturn(0);
1615c66b693SKris Buschelman }
1621c24bd37SHong Zhang 
1631c24bd37SHong Zhang EXTERN int MatDestroy_MPIAIJ(Mat);
16426be0446SHong Zhang #undef __FUNCT__
16526be0446SHong Zhang #define __FUNCT__ "MatDestroy_MPIAIJ_MatMatMult"
16626be0446SHong Zhang int MatDestroy_MPIAIJ_MatMatMult(Mat A)
16726be0446SHong Zhang {
16826be0446SHong Zhang   int               ierr;
16926be0446SHong Zhang   Mat_MatMatMultMPI *mult=(Mat_MatMatMultMPI*)A->spptr;
17026be0446SHong Zhang 
17126be0446SHong Zhang   PetscFunctionBegin;
172d6bb3c2dSHong Zhang   ierr = ISDestroy(mult->isrowb);CHKERRQ(ierr);
173d6bb3c2dSHong Zhang   ierr = ISDestroy(mult->iscolb);CHKERRQ(ierr);
174d6bb3c2dSHong Zhang   ierr = ISDestroy(mult->isrowa);CHKERRQ(ierr);
175d6bb3c2dSHong Zhang   ierr = MatDestroyMatrices(1,&mult->aseq);CHKERRQ(ierr);
176d6bb3c2dSHong Zhang   ierr = MatDestroyMatrices(1,&mult->bseq);CHKERRQ(ierr);
177d6bb3c2dSHong Zhang   ierr = MatDestroy(mult->C_seq);CHKERRQ(ierr);
17826be0446SHong Zhang   ierr = PetscFree(mult);CHKERRQ(ierr);
179d6bb3c2dSHong Zhang 
18026be0446SHong Zhang   ierr = MatDestroy_MPIAIJ(A);CHKERRQ(ierr);
18126be0446SHong Zhang 
18226be0446SHong Zhang   PetscFunctionReturn(0);
18326be0446SHong Zhang }
18458c24d83SHong Zhang 
18558c24d83SHong Zhang #undef __FUNCT__
18626be0446SHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ"
18726be0446SHong Zhang int MatMatMultSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
18826be0446SHong Zhang {
18926be0446SHong Zhang   Mat_MPIAIJ        *a = (Mat_MPIAIJ*)A->data;
19026be0446SHong Zhang   int               ierr,*idx,i,start,end,ncols,imark,nzA,nzB,*cmap;
19126be0446SHong Zhang   Mat_MatMatMultMPI *mult;
19226be0446SHong Zhang 
19326be0446SHong Zhang   PetscFunctionBegin;
194d6bb3c2dSHong Zhang   ierr = PetscNew(Mat_MatMatMultMPI,&mult);CHKERRQ(ierr);
195d6bb3c2dSHong Zhang 
19626be0446SHong Zhang   /* create a seq matrix B_seq = submatrix of B by taking rows of B that equal to nonzero col of A */
19726be0446SHong Zhang   start = a->cstart;
19826be0446SHong Zhang   cmap  = a->garray;
19926be0446SHong Zhang   nzA   = a->A->n;
20026be0446SHong Zhang   nzB   = a->B->n;
20126be0446SHong Zhang   ierr  = PetscMalloc((nzA+nzB)*sizeof(int), &idx);CHKERRQ(ierr);
20226be0446SHong Zhang   ncols = 0;
20326be0446SHong Zhang   for (i=0; i<nzB; i++) {
20426be0446SHong Zhang     if (cmap[i] < start) idx[ncols++] = cmap[i];
20526be0446SHong Zhang     else break;
20626be0446SHong Zhang   }
20726be0446SHong Zhang   imark = i;
20826be0446SHong Zhang   for (i=0; i<nzA; i++) idx[ncols++] = start + i;
20926be0446SHong Zhang   for (i=imark; i<nzB; i++) idx[ncols++] = cmap[i];
210d6bb3c2dSHong Zhang   ierr = ISCreateGeneral(PETSC_COMM_SELF,ncols,idx,&mult->isrowb);CHKERRQ(ierr);
21126be0446SHong Zhang   ierr = PetscFree(idx);CHKERRQ(ierr);
212d6bb3c2dSHong Zhang   ierr = ISCreateStride(PETSC_COMM_SELF,B->N,0,1,&mult->iscolb);CHKERRQ(ierr);
213d6bb3c2dSHong Zhang   ierr = MatGetSubMatrices(B,1,&mult->isrowb,&mult->iscolb,MAT_INITIAL_MATRIX,&mult->bseq);CHKERRQ(ierr)
21426be0446SHong Zhang 
21526be0446SHong Zhang   /*  create a seq matrix A_seq = submatrix of A by taking all local rows of A */
21626be0446SHong Zhang   start = a->rstart; end = a->rend;
217d6bb3c2dSHong Zhang   ierr = ISCreateStride(PETSC_COMM_SELF,end-start,start,1,&mult->isrowa);CHKERRQ(ierr);
218d6bb3c2dSHong Zhang   ierr = MatGetSubMatrices(A,1,&mult->isrowa,&mult->isrowb,MAT_INITIAL_MATRIX,&mult->aseq);CHKERRQ(ierr);
21926be0446SHong Zhang 
22026be0446SHong Zhang   /* compute C_seq = A_seq * B_seq */
221d6bb3c2dSHong Zhang   ierr = MatMatMult_SeqAIJ_SeqAIJ(mult->aseq[0],mult->bseq[0],MAT_INITIAL_MATRIX,fill,&mult->C_seq);CHKERRQ(ierr);
22226be0446SHong Zhang 
22326be0446SHong Zhang   /* create mpi matrix C by concatinating C_seq */
224d6bb3c2dSHong Zhang   ierr = PetscObjectReference((PetscObject)mult->C_seq);CHKERRQ(ierr); /* prevent C_seq being destroyed by MatMerge() */
225d6bb3c2dSHong Zhang   ierr = MatMerge(A->comm,mult->C_seq,MAT_INITIAL_MATRIX,C);CHKERRQ(ierr);
22626be0446SHong Zhang 
22726be0446SHong Zhang   /* attach the supporting struct to C for reuse of symbolic C */
22826be0446SHong Zhang   (*C)->spptr         = (void*)mult;
22926be0446SHong Zhang   (*C)->ops->destroy  = MatDestroy_MPIAIJ_MatMatMult;
23026be0446SHong Zhang 
23126be0446SHong Zhang   PetscFunctionReturn(0);
23226be0446SHong Zhang }
23326be0446SHong Zhang 
23426be0446SHong Zhang #undef __FUNCT__
23526be0446SHong Zhang #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ"
23626be0446SHong Zhang int MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
23758c24d83SHong Zhang {
23858c24d83SHong Zhang   int            ierr;
23958c24d83SHong Zhang   FreeSpaceList  free_space=PETSC_NULL,current_space=PETSC_NULL;
24058c24d83SHong Zhang   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
24158c24d83SHong Zhang   int            *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj;
2421c24bd37SHong Zhang   int            *ci,*cj,*lnk;
2435c66b693SKris Buschelman   int            am=A->M,bn=B->N,bm=B->M;
244c5db241fSHong Zhang   int            i,j,anzi,brow,bnzj,cnzi,nlnk,lnk_init=-1,nspacedouble=0;
24558c24d83SHong Zhang   MatScalar      *ca;
24658c24d83SHong Zhang 
24758c24d83SHong Zhang   PetscFunctionBegin;
2485c66b693SKris Buschelman   /* Start timers */
24958c24d83SHong Zhang   ierr = PetscLogEventBegin(logkey_matmatmult_symbolic,A,B,0,0);CHKERRQ(ierr);
25058c24d83SHong Zhang   /* Set up */
25158c24d83SHong Zhang   /* Allocate ci array, arrays for fill computation and */
25258c24d83SHong Zhang   /* free space for accumulating nonzero column info */
25358c24d83SHong Zhang   ierr = PetscMalloc(((am+1)+1)*sizeof(int),&ci);CHKERRQ(ierr);
25458c24d83SHong Zhang   ci[0] = 0;
25558c24d83SHong Zhang 
25658c24d83SHong Zhang   ierr = PetscMalloc((bn+1)*sizeof(int),&lnk);CHKERRQ(ierr);
257*f747e1aeSHong Zhang   nlnk = bn+1;
258*f747e1aeSHong Zhang   ierr = PetscLLInitialize(lnk_init,nlnk,lnk);CHKERRQ(ierr);
25958c24d83SHong Zhang 
260c5db241fSHong Zhang   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
261d6bb3c2dSHong Zhang   ierr = GetMoreSpace((int)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr);
26258c24d83SHong Zhang   current_space = free_space;
26358c24d83SHong Zhang 
26458c24d83SHong Zhang   /* Determine symbolic info for each row of the product: */
26558c24d83SHong Zhang   for (i=0;i<am;i++) {
26658c24d83SHong Zhang     anzi = ai[i+1] - ai[i];
26758c24d83SHong Zhang     cnzi = 0;
26858c24d83SHong Zhang     lnk[bn] = bn;
2692d09714cSHong Zhang     j       = anzi;
2702d09714cSHong Zhang     aj      = a->j + ai[i];
2712d09714cSHong Zhang     while (j){/* assume cols are almost in increasing order, starting from its end saves computation */
2722d09714cSHong Zhang       j--;
2732d09714cSHong Zhang       brow = *(aj + j);
27458c24d83SHong Zhang       bnzj = bi[brow+1] - bi[brow];
27558c24d83SHong Zhang       bjj  = bj + bi[brow];
2761c239cc6SHong Zhang       /* add non-zero cols of B into the sorted linked list lnk */
27745a8bf62SHong Zhang       ierr = PetscLLAdd(bnzj,bjj,bn,lnk_init,nlnk,lnk);CHKERRQ(ierr);
2781c239cc6SHong Zhang       cnzi += nlnk;
27958c24d83SHong Zhang     }
28058c24d83SHong Zhang 
28158c24d83SHong Zhang     /* If free space is not available, make more free space */
28258c24d83SHong Zhang     /* Double the amount of total space in the list */
28358c24d83SHong Zhang     if (current_space->local_remaining<cnzi) {
28426be0446SHong Zhang       /* printf("MatMatMultSymbolic_SeqAIJ_SeqAIJ()...%d -th row, double space ...\n",i);*/
28558c24d83SHong Zhang       ierr = GetMoreSpace(current_space->total_array_size,&current_space);CHKERRQ(ierr);
286c5db241fSHong Zhang       nspacedouble++;
28758c24d83SHong Zhang     }
28858c24d83SHong Zhang 
289c5db241fSHong Zhang     /* Copy data into free space, then initialize lnk */
29045a8bf62SHong Zhang     ierr = PetscLLClear(bn,lnk_init,cnzi,lnk,current_space->array);CHKERRQ(ierr);
291c5db241fSHong Zhang     current_space->array += cnzi;
29258c24d83SHong Zhang 
29358c24d83SHong Zhang     current_space->local_used      += cnzi;
29458c24d83SHong Zhang     current_space->local_remaining -= cnzi;
29558c24d83SHong Zhang 
29658c24d83SHong Zhang     ci[i+1] = ci[i] + cnzi;
29758c24d83SHong Zhang   }
29858c24d83SHong Zhang 
29958c24d83SHong Zhang   /* Column indices are in the list of free space */
30058c24d83SHong Zhang   /* Allocate space for cj, initialize cj, and */
30158c24d83SHong Zhang   /* destroy list of free space and other temporary array(s) */
30258c24d83SHong Zhang   ierr = PetscMalloc((ci[am]+1)*sizeof(int),&cj);CHKERRQ(ierr);
30358c24d83SHong Zhang   ierr = MakeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
30458c24d83SHong Zhang   ierr = PetscFree(lnk);CHKERRQ(ierr);
30558c24d83SHong Zhang 
30658c24d83SHong Zhang   /* Allocate space for ca */
30758c24d83SHong Zhang   ierr = PetscMalloc((ci[am]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
30858c24d83SHong Zhang   ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr);
30958c24d83SHong Zhang 
31026be0446SHong Zhang   /* put together the new symbolic matrix */
31158c24d83SHong Zhang   ierr = MatCreateSeqAIJWithArrays(A->comm,am,bn,ci,cj,ca,C);CHKERRQ(ierr);
31258c24d83SHong Zhang 
31358c24d83SHong Zhang   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
31458c24d83SHong Zhang   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
31558c24d83SHong Zhang   c = (Mat_SeqAIJ *)((*C)->data);
31658c24d83SHong Zhang   c->freedata = PETSC_TRUE;
31758c24d83SHong Zhang   c->nonew    = 0;
31858c24d83SHong Zhang 
31958c24d83SHong Zhang   ierr = PetscLogEventEnd(logkey_matmatmult_symbolic,A,B,0,0);CHKERRQ(ierr);
320c5db241fSHong Zhang   PetscLogInfo((PetscObject)(*C),"Number of calls to GetMoreSpace(): %d\n",nspacedouble);
32158c24d83SHong Zhang   PetscFunctionReturn(0);
32258c24d83SHong Zhang }
323d50806bdSBarry Smith 
324c1f4806aSKris Buschelman #undef __FUNCT__
325c1f4806aSKris Buschelman #define __FUNCT__ "MatMatMultNumeric"
3265c66b693SKris Buschelman /*@
3275c66b693SKris Buschelman    MatMatMultNumeric - Performs the numeric matrix-matrix product.
3285c66b693SKris Buschelman    Call this routine after first calling MatMatMultSymbolic().
3295c66b693SKris Buschelman 
3305c66b693SKris Buschelman    Collective on Mat
3315c66b693SKris Buschelman 
3325c66b693SKris Buschelman    Input Parameters:
3335c66b693SKris Buschelman +  A - the left matrix
3345c66b693SKris Buschelman -  B - the right matrix
3355c66b693SKris Buschelman 
3365c66b693SKris Buschelman    Output Parameters:
3375c66b693SKris Buschelman .  C - the product matrix, whose ij structure was defined from MatMatMultSymbolic().
3385c66b693SKris Buschelman 
3395c66b693SKris Buschelman    Notes:
3405c66b693SKris Buschelman    C must have been created with MatMatMultSymbolic.
3415c66b693SKris Buschelman 
3425c66b693SKris Buschelman    This routine is currently only implemented for SeqAIJ type matrices.
3435c66b693SKris Buschelman 
3445c66b693SKris Buschelman    Level: intermediate
3455c66b693SKris Buschelman 
3465c66b693SKris Buschelman .seealso: MatMatMult(),MatMatMultSymbolic()
3475c66b693SKris Buschelman @*/
3485c66b693SKris Buschelman int MatMatMultNumeric(Mat A,Mat B,Mat C){
3495c66b693SKris Buschelman   /* Perhaps this "interface" routine should be moved into the interface directory.*/
3505c66b693SKris Buschelman   /* To facilitate implementations with varying types, QueryFunction is used.*/
3515c66b693SKris Buschelman   /* It is assumed that implementations will be composed as "MatMatMultNumeric_<type of A><type of B>". */
3525c66b693SKris Buschelman   int ierr;
3534d3841fdSKris Buschelman   char numfunct[80];
3545c66b693SKris Buschelman   int (*numeric)(Mat,Mat,Mat);
3555c66b693SKris Buschelman 
3565c66b693SKris Buschelman   PetscFunctionBegin;
3575c66b693SKris Buschelman 
3584482741eSBarry Smith   PetscValidHeaderSpecific(A,MAT_COOKIE,1);
359c9780b6fSBarry Smith   PetscValidType(A,1);
3605c66b693SKris Buschelman   MatPreallocated(A);
3615c66b693SKris Buschelman   if (!A->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
3625c66b693SKris Buschelman   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
3635c66b693SKris Buschelman 
3644482741eSBarry Smith   PetscValidHeaderSpecific(B,MAT_COOKIE,2);
365c9780b6fSBarry Smith   PetscValidType(B,2);
3665c66b693SKris Buschelman   MatPreallocated(B);
3675c66b693SKris Buschelman   if (!B->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
3685c66b693SKris Buschelman   if (B->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
3695c66b693SKris Buschelman 
3704482741eSBarry Smith   PetscValidHeaderSpecific(C,MAT_COOKIE,3);
371c9780b6fSBarry Smith   PetscValidType(C,3);
3725c66b693SKris Buschelman   MatPreallocated(C);
3735c66b693SKris Buschelman   if (!C->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
3745c66b693SKris Buschelman   if (C->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
3755c66b693SKris Buschelman 
3765c66b693SKris Buschelman   if (B->N!=C->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",B->N,C->N);
3775c66b693SKris Buschelman   if (B->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",B->M,A->N);
3785c66b693SKris Buschelman   if (A->M!=C->M) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",A->M,C->M);
3795c66b693SKris Buschelman 
3804d3841fdSKris Buschelman   /* Currently only _seqaijseqaij is implemented, so just query for it in A and B. */
3814d3841fdSKris Buschelman   /* When other implementations exist, attack the multiple dispatch problem. */
3824d3841fdSKris Buschelman   ierr = PetscStrcpy(numfunct,"MatMatMultNumeric_seqaijseqaij");CHKERRQ(ierr);
3834d3841fdSKris Buschelman   ierr = PetscObjectQueryFunction((PetscObject)A,numfunct,(PetscVoidFunction)&numeric);CHKERRQ(ierr);
3844d3841fdSKris Buschelman   if (!numeric) SETERRQ1(PETSC_ERR_SUP,"C=A*B not implemented for A of type %s",A->type_name);
3854d3841fdSKris Buschelman   ierr = PetscObjectQueryFunction((PetscObject)B,numfunct,(PetscVoidFunction)&numeric);CHKERRQ(ierr);
3864d3841fdSKris Buschelman   if (!numeric) SETERRQ1(PETSC_ERR_SUP,"C=A*B not implemented for B of type %s",B->type_name);
3874d3841fdSKris Buschelman 
3885c66b693SKris Buschelman   ierr = (*numeric)(A,B,C);CHKERRQ(ierr);
3895c66b693SKris Buschelman 
3905c66b693SKris Buschelman   PetscFunctionReturn(0);
3915c66b693SKris Buschelman }
3925c66b693SKris Buschelman 
393d6bb3c2dSHong Zhang /* This routine is called ONLY in the case of reusing previously computed symbolic C */
394d50806bdSBarry Smith #undef __FUNCT__
39526be0446SHong Zhang #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIAIJ"
39626be0446SHong Zhang int MatMatMultNumeric_MPIAIJ_MPIAIJ(Mat A,Mat B,Mat C)
39726be0446SHong Zhang {
398d6bb3c2dSHong Zhang   int               ierr;
399d6bb3c2dSHong Zhang   Mat_MatMatMultMPI *mult=(Mat_MatMatMultMPI*)C->spptr;
400d6bb3c2dSHong Zhang 
40126be0446SHong Zhang   PetscFunctionBegin;
402d6bb3c2dSHong Zhang   ierr = MatGetSubMatrices(B,1,&mult->isrowb,&mult->iscolb,MAT_REUSE_MATRIX,&mult->bseq);CHKERRQ(ierr)
403d6bb3c2dSHong Zhang   ierr = MatGetSubMatrices(A,1,&mult->isrowa,&mult->isrowb,MAT_REUSE_MATRIX,&mult->aseq);CHKERRQ(ierr);
404d6bb3c2dSHong Zhang   ierr = MatMatMult_SeqAIJ_SeqAIJ(mult->aseq[0],mult->bseq[0],MAT_REUSE_MATRIX,0.0,&mult->C_seq);CHKERRQ(ierr);
405d6bb3c2dSHong Zhang 
406d6bb3c2dSHong Zhang   ierr = PetscObjectReference((PetscObject)mult->C_seq);CHKERRQ(ierr);
407d6bb3c2dSHong Zhang   ierr = MatMerge(A->comm,mult->C_seq,MAT_REUSE_MATRIX,&C);CHKERRQ(ierr);
408d6bb3c2dSHong Zhang 
40926be0446SHong Zhang   PetscFunctionReturn(0);
41026be0446SHong Zhang }
41126be0446SHong Zhang 
41226be0446SHong Zhang #undef __FUNCT__
41326be0446SHong Zhang #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqAIJ"
41426be0446SHong Zhang int MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
415d50806bdSBarry Smith {
41694e3eecaSKris Buschelman   int        ierr,flops=0;
417d50806bdSBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data;
418d50806bdSBarry Smith   Mat_SeqAIJ *b = (Mat_SeqAIJ *)B->data;
419d50806bdSBarry Smith   Mat_SeqAIJ *c = (Mat_SeqAIJ *)C->data;
420d50806bdSBarry Smith   int        *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j;
4215c66b693SKris Buschelman   int        am=A->M,cn=C->N;
42294e3eecaSKris Buschelman   int        i,j,k,anzi,bnzi,cnzi,brow;
423d50806bdSBarry Smith   MatScalar  *aa=a->a,*ba=b->a,*baj,*ca=c->a,*temp;
424d50806bdSBarry Smith 
425d50806bdSBarry Smith   PetscFunctionBegin;
426d50806bdSBarry Smith 
4275c66b693SKris Buschelman   /* Start timers */
428d50806bdSBarry Smith   ierr = PetscLogEventBegin(logkey_matmatmult_numeric,A,B,C,0);CHKERRQ(ierr);
42994e3eecaSKris Buschelman 
430d50806bdSBarry Smith   /* Allocate temp accumulation space to avoid searching for nonzero columns in C */
431d50806bdSBarry Smith   ierr = PetscMalloc((cn+1)*sizeof(MatScalar),&temp);CHKERRQ(ierr);
432d50806bdSBarry Smith   ierr = PetscMemzero(temp,cn*sizeof(MatScalar));CHKERRQ(ierr);
433d50806bdSBarry Smith   /* Traverse A row-wise. */
434d50806bdSBarry Smith   /* Build the ith row in C by summing over nonzero columns in A, */
435d50806bdSBarry Smith   /* the rows of B corresponding to nonzeros of A. */
436d50806bdSBarry Smith   for (i=0;i<am;i++) {
437d50806bdSBarry Smith     anzi = ai[i+1] - ai[i];
438d50806bdSBarry Smith     for (j=0;j<anzi;j++) {
439d50806bdSBarry Smith       brow = *aj++;
440d50806bdSBarry Smith       bnzi = bi[brow+1] - bi[brow];
441d50806bdSBarry Smith       bjj  = bj + bi[brow];
442d50806bdSBarry Smith       baj  = ba + bi[brow];
443d50806bdSBarry Smith       for (k=0;k<bnzi;k++) {
444d50806bdSBarry Smith         temp[bjj[k]] += (*aa)*baj[k];
445d50806bdSBarry Smith       }
446d50806bdSBarry Smith       flops += 2*bnzi;
447d50806bdSBarry Smith       aa++;
448d50806bdSBarry Smith     }
449d50806bdSBarry Smith     /* Store row back into C, and re-zero temp */
450d50806bdSBarry Smith     cnzi = ci[i+1] - ci[i];
451d50806bdSBarry Smith     for (j=0;j<cnzi;j++) {
452d50806bdSBarry Smith       ca[j] = temp[cj[j]];
453d50806bdSBarry Smith       temp[cj[j]] = 0.0;
454d50806bdSBarry Smith     }
455d50806bdSBarry Smith     ca += cnzi;
456d50806bdSBarry Smith     cj += cnzi;
457d50806bdSBarry Smith   }
458716bacf3SKris Buschelman   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
459716bacf3SKris Buschelman   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
460716bacf3SKris Buschelman 
461d50806bdSBarry Smith   /* Free temp */
462d50806bdSBarry Smith   ierr = PetscFree(temp);CHKERRQ(ierr);
463d50806bdSBarry Smith   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
464d50806bdSBarry Smith   ierr = PetscLogEventEnd(logkey_matmatmult_numeric,A,B,C,0);CHKERRQ(ierr);
465d50806bdSBarry Smith   PetscFunctionReturn(0);
466d50806bdSBarry Smith }
467d50806bdSBarry Smith 
468d50806bdSBarry Smith #undef __FUNCT__
4695c66b693SKris Buschelman #define __FUNCT__ "RegisterMatMatMultRoutines_Private"
4705c66b693SKris Buschelman int RegisterMatMatMultRoutines_Private(Mat A) {
471d50806bdSBarry Smith   int ierr;
472d50806bdSBarry Smith 
473d50806bdSBarry Smith   PetscFunctionBegin;
4745c66b693SKris Buschelman   if (!logkey_matmatmult_symbolic) {
47526be0446SHong Zhang     ierr = PetscLogEventRegister(&logkey_matmatmult_symbolic,"MatMatMultSymbolic",MAT_COOKIE);CHKERRQ(ierr);
476d50806bdSBarry Smith   }
4775c66b693SKris Buschelman   if (!logkey_matmatmult_numeric) {
47826be0446SHong Zhang     ierr = PetscLogEventRegister(&logkey_matmatmult_numeric,"MatMatMultNumeric",MAT_COOKIE);CHKERRQ(ierr);
47994e3eecaSKris Buschelman   }
480d6bb3c2dSHong Zhang 
48194e3eecaSKris Buschelman   PetscFunctionReturn(0);
48294e3eecaSKris Buschelman }
483