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 15c1f4806aSKris Buschelman #undef __FUNCT__ 16c1f4806aSKris Buschelman #define __FUNCT__ "MatMatMult" 175c66b693SKris Buschelman /*@ 185c66b693SKris Buschelman MatMatMult - Performs Matrix-Matrix Multiplication C=A*B. 1994e3eecaSKris Buschelman 205c66b693SKris Buschelman Collective on Mat 21d50806bdSBarry Smith 225c66b693SKris Buschelman Input Parameters: 235c66b693SKris Buschelman + A - the left matrix 241c741599SHong Zhang . B - the right matrix 251c741599SHong Zhang . scall - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX 26c5db241fSHong Zhang - fill - expected fill as ratio of nnz(C)/(nnz(A) + nnz(B)) 275c66b693SKris Buschelman 285c66b693SKris Buschelman Output Parameters: 295c66b693SKris Buschelman . C - the product matrix 305c66b693SKris Buschelman 315c66b693SKris Buschelman Notes: 325c66b693SKris Buschelman C will be created and must be destroyed by the user with MatDestroy(). 335c66b693SKris Buschelman 344d3841fdSKris Buschelman This routine is currently only implemented for pairs of SeqAIJ matrices and classes 354d3841fdSKris Buschelman which inherit from SeqAIJ. C will be of type MATSEQAIJ. 365c66b693SKris Buschelman 375c66b693SKris Buschelman Level: intermediate 385c66b693SKris Buschelman 395c66b693SKris Buschelman .seealso: MatMatMultSymbolic(),MatMatMultNumeric() 405c66b693SKris Buschelman @*/ 41dfbe8321SBarry Smith PetscErrorCode MatMatMult(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 422d09714cSHong Zhang { 43dfbe8321SBarry Smith PetscErrorCode ierr; 44cb00f407SKris Buschelman PetscErrorCode (*fA)(Mat,Mat,MatReuse,PetscReal,Mat*); 45cb00f407SKris Buschelman PetscErrorCode (*fB)(Mat,Mat,MatReuse,PetscReal,Mat*); 462d09714cSHong Zhang 47d50806bdSBarry Smith PetscFunctionBegin; 484482741eSBarry Smith PetscValidHeaderSpecific(A,MAT_COOKIE,1); 49c9780b6fSBarry Smith PetscValidType(A,1); 505c66b693SKris Buschelman MatPreallocated(A); 515c66b693SKris Buschelman if (!A->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); 525c66b693SKris Buschelman if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 534482741eSBarry Smith PetscValidHeaderSpecific(B,MAT_COOKIE,2); 54c9780b6fSBarry Smith PetscValidType(B,2); 555c66b693SKris Buschelman MatPreallocated(B); 565c66b693SKris Buschelman if (!B->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); 575c66b693SKris Buschelman if (B->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 584482741eSBarry Smith PetscValidPointer(C,3); 595c66b693SKris Buschelman if (B->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",B->M,A->N); 60d50806bdSBarry Smith 611c24bd37SHong Zhang if (fill <=0.0) SETERRQ1(PETSC_ERR_ARG_SIZ,"fill=%g must be > 0.0",fill); 62c5db241fSHong Zhang 63cb00f407SKris Buschelman /* For now, we do not dispatch based on the type of A and B */ 64cb00f407SKris Buschelman /* When implementations like _SeqAIJ_MAIJ exist, attack the multiple dispatch problem. */ 65cb00f407SKris Buschelman fA = A->ops->matmult; 66cb00f407SKris Buschelman if (!fA) SETERRQ1(PETSC_ERR_SUP,"MatMatMult not supported for A of type %s",A->type_name); 67cb00f407SKris Buschelman fB = B->ops->matmult; 68cb00f407SKris Buschelman if (!fB) SETERRQ1(PETSC_ERR_SUP,"MatMatMult not supported for B of type %s",B->type_name); 69cb00f407SKris Buschelman if (fB!=fA) SETERRQ2(PETSC_ERR_ARG_INCOMP,"MatMatMult requires A, %s, to be compatible with B, %s",A->type_name,B->type_name); 70cb00f407SKris Buschelman 716284ec50SHong Zhang ierr = PetscLogEventBegin(MAT_MatMult,A,B,0,0);CHKERRQ(ierr); 721c741599SHong Zhang ierr = (*A->ops->matmult)(A,B,scall,fill,C);CHKERRQ(ierr); 736284ec50SHong Zhang ierr = PetscLogEventEnd(MAT_MatMult,A,B,0,0);CHKERRQ(ierr); 744d3841fdSKris Buschelman 756284ec50SHong Zhang PetscFunctionReturn(0); 766284ec50SHong Zhang } 776284ec50SHong Zhang 786284ec50SHong Zhang #undef __FUNCT__ 796284ec50SHong Zhang #define __FUNCT__ "MatMatMult_MPIAIJ_MPIAIJ" 80dfbe8321SBarry Smith PetscErrorCode MatMatMult_MPIAIJ_MPIAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill, Mat *C) 812d09714cSHong Zhang { 82dfbe8321SBarry Smith PetscErrorCode ierr; 836284ec50SHong Zhang 846284ec50SHong Zhang PetscFunctionBegin; 8526be0446SHong Zhang if (scall == MAT_INITIAL_MATRIX){ 86d6bb3c2dSHong Zhang ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ(A,B,fill,C);CHKERRQ(ierr);/* numeric product is computed as well */ 87d6bb3c2dSHong Zhang } else if (scall == MAT_REUSE_MATRIX){ 8826be0446SHong Zhang ierr = MatMatMultNumeric_MPIAIJ_MPIAIJ(A,B,*C);CHKERRQ(ierr); 89d6bb3c2dSHong Zhang } else { 90d6bb3c2dSHong Zhang SETERRQ1(PETSC_ERR_ARG_WRONG,"Invalid MatReuse %d",scall); 91d6bb3c2dSHong Zhang } 92d50806bdSBarry Smith PetscFunctionReturn(0); 93d50806bdSBarry Smith } 945c66b693SKris Buschelman 955c66b693SKris Buschelman #undef __FUNCT__ 965c66b693SKris Buschelman #define __FUNCT__ "MatMatMult_SeqAIJ_SeqAIJ" 97dfbe8321SBarry Smith PetscErrorCode MatMatMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) { 98dfbe8321SBarry Smith PetscErrorCode ierr; 995c66b693SKris Buschelman 1005c66b693SKris Buschelman PetscFunctionBegin; 10126be0446SHong Zhang if (scall == MAT_INITIAL_MATRIX){ 10226be0446SHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr); 10326be0446SHong Zhang } 10426be0446SHong Zhang ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr); 1055c66b693SKris Buschelman PetscFunctionReturn(0); 1065c66b693SKris Buschelman } 1075c66b693SKris Buschelman 108c1f4806aSKris Buschelman #undef __FUNCT__ 109c1f4806aSKris Buschelman #define __FUNCT__ "MatMatMultSymbolic" 1105c66b693SKris Buschelman /*@ 1115c66b693SKris Buschelman MatMatMultSymbolic - Performs construction, preallocation, and computes the ij structure 1125c66b693SKris Buschelman of the matrix-matrix product C=A*B. Call this routine before calling MatMatMultNumeric(). 1135c66b693SKris Buschelman 1145c66b693SKris Buschelman Collective on Mat 1155c66b693SKris Buschelman 1165c66b693SKris Buschelman Input Parameters: 1175c66b693SKris Buschelman + A - the left matrix 118c5db241fSHong Zhang . B - the right matrix 119c5db241fSHong Zhang - fill - expected fill as ratio of nnz(C)/(nnz(A) + nnz(B)) 1205c66b693SKris Buschelman 1215c66b693SKris Buschelman Output Parameters: 1225c66b693SKris Buschelman . C - the matrix containing the ij structure of product matrix 1235c66b693SKris Buschelman 1245c66b693SKris Buschelman Notes: 1254d3841fdSKris Buschelman C will be created as a MATSEQAIJ matrix and must be destroyed by the user with MatDestroy(). 1265c66b693SKris Buschelman 1274d3841fdSKris Buschelman This routine is currently only implemented for SeqAIJ matrices and classes which inherit from SeqAIJ. 1285c66b693SKris Buschelman 1295c66b693SKris Buschelman Level: intermediate 1305c66b693SKris Buschelman 1315c66b693SKris Buschelman .seealso: MatMatMult(),MatMatMultNumeric() 1325c66b693SKris Buschelman @*/ 133dfbe8321SBarry Smith PetscErrorCode MatMatMultSymbolic(Mat A,Mat B,PetscReal fill,Mat *C) { 134dfbe8321SBarry Smith PetscErrorCode ierr; 135cb00f407SKris Buschelman PetscErrorCode (*Asymbolic)(Mat,Mat,PetscReal,Mat *); 136cb00f407SKris Buschelman PetscErrorCode (*Bsymbolic)(Mat,Mat,PetscReal,Mat *); 1375c66b693SKris Buschelman 1385c66b693SKris Buschelman PetscFunctionBegin; 1394482741eSBarry Smith PetscValidHeaderSpecific(A,MAT_COOKIE,1); 140c9780b6fSBarry Smith PetscValidType(A,1); 1415c66b693SKris Buschelman MatPreallocated(A); 1425c66b693SKris Buschelman if (!A->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); 1435c66b693SKris Buschelman if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 1445c66b693SKris Buschelman 1454482741eSBarry Smith PetscValidHeaderSpecific(B,MAT_COOKIE,2); 146c9780b6fSBarry Smith PetscValidType(B,2); 1475c66b693SKris Buschelman MatPreallocated(B); 1485c66b693SKris Buschelman if (!B->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); 1495c66b693SKris Buschelman if (B->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 1504482741eSBarry Smith PetscValidPointer(C,3); 1514482741eSBarry Smith 1525c66b693SKris Buschelman if (B->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",B->M,A->N); 1531c24bd37SHong Zhang if (fill <=0.0) SETERRQ1(PETSC_ERR_ARG_SIZ,"fill=%g must be > 0.0",fill); 1545c66b693SKris Buschelman 155cb00f407SKris Buschelman /* For now, we do not dispatch based on the type of A and P */ 156cb00f407SKris Buschelman /* When implementations like _SeqAIJ_MAIJ exist, attack the multiple dispatch problem. */ 157cb00f407SKris Buschelman Asymbolic = A->ops->matmultsymbolic; 158cb00f407SKris Buschelman if (!Asymbolic) SETERRQ1(PETSC_ERR_SUP,"C=A*B not implemented for A of type %s",A->type_name); 159cb00f407SKris Buschelman Bsymbolic = B->ops->matmultsymbolic; 160cb00f407SKris Buschelman if (!Bsymbolic) SETERRQ1(PETSC_ERR_SUP,"C=A*B not implemented for B of type %s",B->type_name); 161cb00f407SKris Buschelman if (Bsymbolic!=Asymbolic) SETERRQ2(PETSC_ERR_ARG_INCOMP,"MatMatMultSymbolic requires A, %s, to be compatible with B, %s",A->type_name,B->type_name); 162cb00f407SKris Buschelman 163cb00f407SKris Buschelman ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr); 164cb00f407SKris Buschelman ierr = (*Asymbolic)(A,B,fill,C);CHKERRQ(ierr); 165cb00f407SKris Buschelman ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr); 1665c66b693SKris Buschelman 1675c66b693SKris Buschelman PetscFunctionReturn(0); 1685c66b693SKris Buschelman } 1691c24bd37SHong Zhang 170dfbe8321SBarry Smith EXTERN PetscErrorCode MatDestroy_MPIAIJ(Mat); 17126be0446SHong Zhang #undef __FUNCT__ 17226be0446SHong Zhang #define __FUNCT__ "MatDestroy_MPIAIJ_MatMatMult" 173dfbe8321SBarry Smith PetscErrorCode MatDestroy_MPIAIJ_MatMatMult(Mat A) 17426be0446SHong Zhang { 175dfbe8321SBarry Smith PetscErrorCode ierr; 17626be0446SHong Zhang Mat_MatMatMultMPI *mult=(Mat_MatMatMultMPI*)A->spptr; 17726be0446SHong Zhang 17826be0446SHong Zhang PetscFunctionBegin; 179d6bb3c2dSHong Zhang ierr = ISDestroy(mult->isrowb);CHKERRQ(ierr); 180d6bb3c2dSHong Zhang ierr = ISDestroy(mult->iscolb);CHKERRQ(ierr); 181d6bb3c2dSHong Zhang ierr = ISDestroy(mult->isrowa);CHKERRQ(ierr); 182d6bb3c2dSHong Zhang ierr = MatDestroyMatrices(1,&mult->aseq);CHKERRQ(ierr); 183d6bb3c2dSHong Zhang ierr = MatDestroyMatrices(1,&mult->bseq);CHKERRQ(ierr); 184d6bb3c2dSHong Zhang ierr = MatDestroy(mult->C_seq);CHKERRQ(ierr); 18526be0446SHong Zhang ierr = PetscFree(mult);CHKERRQ(ierr); 186d6bb3c2dSHong Zhang 18726be0446SHong Zhang ierr = MatDestroy_MPIAIJ(A);CHKERRQ(ierr); 18826be0446SHong Zhang 18926be0446SHong Zhang PetscFunctionReturn(0); 19026be0446SHong Zhang } 19158c24d83SHong Zhang 19258c24d83SHong Zhang #undef __FUNCT__ 19326be0446SHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ" 194dfbe8321SBarry Smith PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat B,PetscReal fill,Mat *C) 19526be0446SHong Zhang { 19626be0446SHong Zhang Mat_MPIAIJ *a = (Mat_MPIAIJ*)A->data; 197dfbe8321SBarry Smith PetscErrorCode ierr; 198dfbe8321SBarry Smith int *idx,i,start,end,ncols,imark,nzA,nzB,*cmap; 19926be0446SHong Zhang Mat_MatMatMultMPI *mult; 20026be0446SHong Zhang 20126be0446SHong Zhang PetscFunctionBegin; 202d6bb3c2dSHong Zhang ierr = PetscNew(Mat_MatMatMultMPI,&mult);CHKERRQ(ierr); 203d6bb3c2dSHong Zhang 20426be0446SHong Zhang /* create a seq matrix B_seq = submatrix of B by taking rows of B that equal to nonzero col of A */ 20526be0446SHong Zhang start = a->cstart; 20626be0446SHong Zhang cmap = a->garray; 20726be0446SHong Zhang nzA = a->A->n; 20826be0446SHong Zhang nzB = a->B->n; 20926be0446SHong Zhang ierr = PetscMalloc((nzA+nzB)*sizeof(int), &idx);CHKERRQ(ierr); 21026be0446SHong Zhang ncols = 0; 21126be0446SHong Zhang for (i=0; i<nzB; i++) { 21226be0446SHong Zhang if (cmap[i] < start) idx[ncols++] = cmap[i]; 21326be0446SHong Zhang else break; 21426be0446SHong Zhang } 21526be0446SHong Zhang imark = i; 21626be0446SHong Zhang for (i=0; i<nzA; i++) idx[ncols++] = start + i; 21726be0446SHong Zhang for (i=imark; i<nzB; i++) idx[ncols++] = cmap[i]; 218d6bb3c2dSHong Zhang ierr = ISCreateGeneral(PETSC_COMM_SELF,ncols,idx,&mult->isrowb);CHKERRQ(ierr); 21926be0446SHong Zhang ierr = PetscFree(idx);CHKERRQ(ierr); 220d6bb3c2dSHong Zhang ierr = ISCreateStride(PETSC_COMM_SELF,B->N,0,1,&mult->iscolb);CHKERRQ(ierr); 221d6bb3c2dSHong Zhang ierr = MatGetSubMatrices(B,1,&mult->isrowb,&mult->iscolb,MAT_INITIAL_MATRIX,&mult->bseq);CHKERRQ(ierr) 22226be0446SHong Zhang 22326be0446SHong Zhang /* create a seq matrix A_seq = submatrix of A by taking all local rows of A */ 22426be0446SHong Zhang start = a->rstart; end = a->rend; 225d6bb3c2dSHong Zhang ierr = ISCreateStride(PETSC_COMM_SELF,end-start,start,1,&mult->isrowa);CHKERRQ(ierr); 226d6bb3c2dSHong Zhang ierr = MatGetSubMatrices(A,1,&mult->isrowa,&mult->isrowb,MAT_INITIAL_MATRIX,&mult->aseq);CHKERRQ(ierr); 22726be0446SHong Zhang 22826be0446SHong Zhang /* compute C_seq = A_seq * B_seq */ 229d6bb3c2dSHong Zhang ierr = MatMatMult_SeqAIJ_SeqAIJ(mult->aseq[0],mult->bseq[0],MAT_INITIAL_MATRIX,fill,&mult->C_seq);CHKERRQ(ierr); 23026be0446SHong Zhang 23126be0446SHong Zhang /* create mpi matrix C by concatinating C_seq */ 232d6bb3c2dSHong Zhang ierr = PetscObjectReference((PetscObject)mult->C_seq);CHKERRQ(ierr); /* prevent C_seq being destroyed by MatMerge() */ 233d6bb3c2dSHong Zhang ierr = MatMerge(A->comm,mult->C_seq,MAT_INITIAL_MATRIX,C);CHKERRQ(ierr); 23426be0446SHong Zhang 23526be0446SHong Zhang /* attach the supporting struct to C for reuse of symbolic C */ 23626be0446SHong Zhang (*C)->spptr = (void*)mult; 23726be0446SHong Zhang (*C)->ops->destroy = MatDestroy_MPIAIJ_MatMatMult; 23826be0446SHong Zhang 23926be0446SHong Zhang PetscFunctionReturn(0); 24026be0446SHong Zhang } 24126be0446SHong Zhang 24226be0446SHong Zhang #undef __FUNCT__ 24326be0446SHong Zhang #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ" 244dfbe8321SBarry Smith PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C) 24558c24d83SHong Zhang { 246dfbe8321SBarry Smith PetscErrorCode ierr; 24758c24d83SHong Zhang FreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 24858c24d83SHong Zhang Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c; 24958c24d83SHong Zhang int *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj; 2501c24bd37SHong Zhang int *ci,*cj,*lnk; 2515c66b693SKris Buschelman int am=A->M,bn=B->N,bm=B->M; 252*0e76890bSHong Zhang int i,j,anzi,brow,bnzj,cnzi,nlnk,nspacedouble=0; 25358c24d83SHong Zhang MatScalar *ca; 25458c24d83SHong Zhang 25558c24d83SHong Zhang PetscFunctionBegin; 25658c24d83SHong Zhang /* Set up */ 25758c24d83SHong Zhang /* Allocate ci array, arrays for fill computation and */ 25858c24d83SHong Zhang /* free space for accumulating nonzero column info */ 25958c24d83SHong Zhang ierr = PetscMalloc(((am+1)+1)*sizeof(int),&ci);CHKERRQ(ierr); 26058c24d83SHong Zhang ci[0] = 0; 26158c24d83SHong Zhang 262*0e76890bSHong Zhang /* create and initialize a linked list for symbolic product */ 26358c24d83SHong Zhang ierr = PetscMalloc((bn+1)*sizeof(int),&lnk);CHKERRQ(ierr); 264*0e76890bSHong Zhang ierr = PetscLLInitialize(bn,bn);CHKERRQ(ierr); 26558c24d83SHong Zhang 266c5db241fSHong Zhang /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */ 267d6bb3c2dSHong Zhang ierr = GetMoreSpace((int)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr); 26858c24d83SHong Zhang current_space = free_space; 26958c24d83SHong Zhang 27058c24d83SHong Zhang /* Determine symbolic info for each row of the product: */ 27158c24d83SHong Zhang for (i=0;i<am;i++) { 27258c24d83SHong Zhang anzi = ai[i+1] - ai[i]; 27358c24d83SHong Zhang cnzi = 0; 2742d09714cSHong Zhang j = anzi; 2752d09714cSHong Zhang aj = a->j + ai[i]; 2762d09714cSHong Zhang while (j){/* assume cols are almost in increasing order, starting from its end saves computation */ 2772d09714cSHong Zhang j--; 2782d09714cSHong Zhang brow = *(aj + j); 27958c24d83SHong Zhang bnzj = bi[brow+1] - bi[brow]; 28058c24d83SHong Zhang bjj = bj + bi[brow]; 2811c239cc6SHong Zhang /* add non-zero cols of B into the sorted linked list lnk */ 282*0e76890bSHong Zhang ierr = PetscLLAdd(bnzj,bjj,bn,nlnk,lnk);CHKERRQ(ierr); 2831c239cc6SHong Zhang cnzi += nlnk; 28458c24d83SHong Zhang } 28558c24d83SHong Zhang 28658c24d83SHong Zhang /* If free space is not available, make more free space */ 28758c24d83SHong Zhang /* Double the amount of total space in the list */ 28858c24d83SHong Zhang if (current_space->local_remaining<cnzi) { 28958c24d83SHong Zhang ierr = GetMoreSpace(current_space->total_array_size,¤t_space);CHKERRQ(ierr); 290c5db241fSHong Zhang nspacedouble++; 29158c24d83SHong Zhang } 29258c24d83SHong Zhang 293c5db241fSHong Zhang /* Copy data into free space, then initialize lnk */ 294*0e76890bSHong Zhang ierr = PetscLLClean(bn,bn,cnzi,lnk,current_space->array);CHKERRQ(ierr); 295c5db241fSHong Zhang current_space->array += cnzi; 29658c24d83SHong Zhang current_space->local_used += cnzi; 29758c24d83SHong Zhang current_space->local_remaining -= cnzi; 29858c24d83SHong Zhang 29958c24d83SHong Zhang ci[i+1] = ci[i] + cnzi; 30058c24d83SHong Zhang } 30158c24d83SHong Zhang 30258c24d83SHong Zhang /* Column indices are in the list of free space */ 30358c24d83SHong Zhang /* Allocate space for cj, initialize cj, and */ 30458c24d83SHong Zhang /* destroy list of free space and other temporary array(s) */ 30558c24d83SHong Zhang ierr = PetscMalloc((ci[am]+1)*sizeof(int),&cj);CHKERRQ(ierr); 30658c24d83SHong Zhang ierr = MakeSpaceContiguous(&free_space,cj);CHKERRQ(ierr); 30758c24d83SHong Zhang ierr = PetscFree(lnk);CHKERRQ(ierr); 30858c24d83SHong Zhang 30958c24d83SHong Zhang /* Allocate space for ca */ 31058c24d83SHong Zhang ierr = PetscMalloc((ci[am]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr); 31158c24d83SHong Zhang ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr); 31258c24d83SHong Zhang 31326be0446SHong Zhang /* put together the new symbolic matrix */ 31458c24d83SHong Zhang ierr = MatCreateSeqAIJWithArrays(A->comm,am,bn,ci,cj,ca,C);CHKERRQ(ierr); 31558c24d83SHong Zhang 31658c24d83SHong Zhang /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 31758c24d83SHong Zhang /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */ 31858c24d83SHong Zhang c = (Mat_SeqAIJ *)((*C)->data); 31958c24d83SHong Zhang c->freedata = PETSC_TRUE; 32058c24d83SHong Zhang c->nonew = 0; 32158c24d83SHong Zhang 322c5db241fSHong Zhang PetscLogInfo((PetscObject)(*C),"Number of calls to GetMoreSpace(): %d\n",nspacedouble); 32358c24d83SHong Zhang PetscFunctionReturn(0); 32458c24d83SHong Zhang } 325d50806bdSBarry Smith 326c1f4806aSKris Buschelman #undef __FUNCT__ 327c1f4806aSKris Buschelman #define __FUNCT__ "MatMatMultNumeric" 3285c66b693SKris Buschelman /*@ 3295c66b693SKris Buschelman MatMatMultNumeric - Performs the numeric matrix-matrix product. 3305c66b693SKris Buschelman Call this routine after first calling MatMatMultSymbolic(). 3315c66b693SKris Buschelman 3325c66b693SKris Buschelman Collective on Mat 3335c66b693SKris Buschelman 3345c66b693SKris Buschelman Input Parameters: 3355c66b693SKris Buschelman + A - the left matrix 3365c66b693SKris Buschelman - B - the right matrix 3375c66b693SKris Buschelman 3385c66b693SKris Buschelman Output Parameters: 3395c66b693SKris Buschelman . C - the product matrix, whose ij structure was defined from MatMatMultSymbolic(). 3405c66b693SKris Buschelman 3415c66b693SKris Buschelman Notes: 3425c66b693SKris Buschelman C must have been created with MatMatMultSymbolic. 3435c66b693SKris Buschelman 3445c66b693SKris Buschelman This routine is currently only implemented for SeqAIJ type matrices. 3455c66b693SKris Buschelman 3465c66b693SKris Buschelman Level: intermediate 3475c66b693SKris Buschelman 3485c66b693SKris Buschelman .seealso: MatMatMult(),MatMatMultSymbolic() 3495c66b693SKris Buschelman @*/ 350dfbe8321SBarry Smith PetscErrorCode MatMatMultNumeric(Mat A,Mat B,Mat C){ 351dfbe8321SBarry Smith PetscErrorCode ierr; 352cb00f407SKris Buschelman PetscErrorCode (*Anumeric)(Mat,Mat,Mat); 353cb00f407SKris Buschelman PetscErrorCode (*Bnumeric)(Mat,Mat,Mat); 3545c66b693SKris Buschelman 3555c66b693SKris Buschelman PetscFunctionBegin; 3565c66b693SKris Buschelman 3574482741eSBarry Smith PetscValidHeaderSpecific(A,MAT_COOKIE,1); 358c9780b6fSBarry Smith PetscValidType(A,1); 3595c66b693SKris Buschelman MatPreallocated(A); 3605c66b693SKris Buschelman if (!A->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); 3615c66b693SKris Buschelman if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 3625c66b693SKris Buschelman 3634482741eSBarry Smith PetscValidHeaderSpecific(B,MAT_COOKIE,2); 364c9780b6fSBarry Smith PetscValidType(B,2); 3655c66b693SKris Buschelman MatPreallocated(B); 3665c66b693SKris Buschelman if (!B->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); 3675c66b693SKris Buschelman if (B->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 3685c66b693SKris Buschelman 3694482741eSBarry Smith PetscValidHeaderSpecific(C,MAT_COOKIE,3); 370c9780b6fSBarry Smith PetscValidType(C,3); 3715c66b693SKris Buschelman MatPreallocated(C); 3725c66b693SKris Buschelman if (!C->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); 3735c66b693SKris Buschelman if (C->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 3745c66b693SKris Buschelman 3755c66b693SKris Buschelman if (B->N!=C->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",B->N,C->N); 3765c66b693SKris Buschelman if (B->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",B->M,A->N); 3775c66b693SKris Buschelman if (A->M!=C->M) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",A->M,C->M); 3785c66b693SKris Buschelman 379cb00f407SKris Buschelman /* For now, we do not dispatch based on the type of A and B */ 380cb00f407SKris Buschelman /* When implementations like _SeqAIJ_MAIJ exist, attack the multiple dispatch problem. */ 381cb00f407SKris Buschelman Anumeric = A->ops->matmultnumeric; 382cb00f407SKris Buschelman if (!Anumeric) SETERRQ1(PETSC_ERR_SUP,"MatMatMultNumeric not supported for A of type %s",A->type_name); 383cb00f407SKris Buschelman Bnumeric = B->ops->matmultnumeric; 384cb00f407SKris Buschelman if (!Bnumeric) SETERRQ1(PETSC_ERR_SUP,"MatMatMultNumeric not supported for B of type %s",B->type_name); 385cb00f407SKris Buschelman if (Bnumeric!=Anumeric) SETERRQ2(PETSC_ERR_ARG_INCOMP,"MatMatMultNumeric requires A, %s, to be compatible with B, %s",A->type_name,B->type_name); 3864d3841fdSKris Buschelman 387cb00f407SKris Buschelman ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr); 388cb00f407SKris Buschelman ierr = (*Anumeric)(A,B,C);CHKERRQ(ierr); 389cb00f407SKris Buschelman ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr); 3905c66b693SKris Buschelman 3915c66b693SKris Buschelman PetscFunctionReturn(0); 3925c66b693SKris Buschelman } 3935c66b693SKris Buschelman 394d6bb3c2dSHong Zhang /* This routine is called ONLY in the case of reusing previously computed symbolic C */ 395d50806bdSBarry Smith #undef __FUNCT__ 39626be0446SHong Zhang #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIAIJ" 397dfbe8321SBarry Smith PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ(Mat A,Mat B,Mat C) 39826be0446SHong Zhang { 399dfbe8321SBarry Smith PetscErrorCode ierr; 400d6bb3c2dSHong Zhang Mat_MatMatMultMPI *mult=(Mat_MatMatMultMPI*)C->spptr; 401d6bb3c2dSHong Zhang 40226be0446SHong Zhang PetscFunctionBegin; 403d6bb3c2dSHong Zhang ierr = MatGetSubMatrices(B,1,&mult->isrowb,&mult->iscolb,MAT_REUSE_MATRIX,&mult->bseq);CHKERRQ(ierr) 404d6bb3c2dSHong Zhang ierr = MatGetSubMatrices(A,1,&mult->isrowa,&mult->isrowb,MAT_REUSE_MATRIX,&mult->aseq);CHKERRQ(ierr); 405d6bb3c2dSHong Zhang ierr = MatMatMult_SeqAIJ_SeqAIJ(mult->aseq[0],mult->bseq[0],MAT_REUSE_MATRIX,0.0,&mult->C_seq);CHKERRQ(ierr); 406d6bb3c2dSHong Zhang 407d6bb3c2dSHong Zhang ierr = PetscObjectReference((PetscObject)mult->C_seq);CHKERRQ(ierr); 408d6bb3c2dSHong Zhang ierr = MatMerge(A->comm,mult->C_seq,MAT_REUSE_MATRIX,&C);CHKERRQ(ierr); 409d6bb3c2dSHong Zhang 41026be0446SHong Zhang PetscFunctionReturn(0); 41126be0446SHong Zhang } 41226be0446SHong Zhang 41326be0446SHong Zhang #undef __FUNCT__ 41426be0446SHong Zhang #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqAIJ" 415dfbe8321SBarry Smith PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C) 416d50806bdSBarry Smith { 417dfbe8321SBarry Smith PetscErrorCode ierr; 418dfbe8321SBarry Smith int flops=0; 419d50806bdSBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 420d50806bdSBarry Smith Mat_SeqAIJ *b = (Mat_SeqAIJ *)B->data; 421d50806bdSBarry Smith Mat_SeqAIJ *c = (Mat_SeqAIJ *)C->data; 422d50806bdSBarry Smith int *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j; 4235c66b693SKris Buschelman int am=A->M,cn=C->N; 42494e3eecaSKris Buschelman int i,j,k,anzi,bnzi,cnzi,brow; 425d50806bdSBarry Smith MatScalar *aa=a->a,*ba=b->a,*baj,*ca=c->a,*temp; 426d50806bdSBarry Smith 427d50806bdSBarry Smith PetscFunctionBegin; 428d50806bdSBarry Smith 429d50806bdSBarry Smith /* Allocate temp accumulation space to avoid searching for nonzero columns in C */ 430d50806bdSBarry Smith ierr = PetscMalloc((cn+1)*sizeof(MatScalar),&temp);CHKERRQ(ierr); 431d50806bdSBarry Smith ierr = PetscMemzero(temp,cn*sizeof(MatScalar));CHKERRQ(ierr); 432d50806bdSBarry Smith /* Traverse A row-wise. */ 433d50806bdSBarry Smith /* Build the ith row in C by summing over nonzero columns in A, */ 434d50806bdSBarry Smith /* the rows of B corresponding to nonzeros of A. */ 435d50806bdSBarry Smith for (i=0;i<am;i++) { 436d50806bdSBarry Smith anzi = ai[i+1] - ai[i]; 437d50806bdSBarry Smith for (j=0;j<anzi;j++) { 438d50806bdSBarry Smith brow = *aj++; 439d50806bdSBarry Smith bnzi = bi[brow+1] - bi[brow]; 440d50806bdSBarry Smith bjj = bj + bi[brow]; 441d50806bdSBarry Smith baj = ba + bi[brow]; 442d50806bdSBarry Smith for (k=0;k<bnzi;k++) { 443d50806bdSBarry Smith temp[bjj[k]] += (*aa)*baj[k]; 444d50806bdSBarry Smith } 445d50806bdSBarry Smith flops += 2*bnzi; 446d50806bdSBarry Smith aa++; 447d50806bdSBarry Smith } 448d50806bdSBarry Smith /* Store row back into C, and re-zero temp */ 449d50806bdSBarry Smith cnzi = ci[i+1] - ci[i]; 450d50806bdSBarry Smith for (j=0;j<cnzi;j++) { 451d50806bdSBarry Smith ca[j] = temp[cj[j]]; 452d50806bdSBarry Smith temp[cj[j]] = 0.0; 453d50806bdSBarry Smith } 454d50806bdSBarry Smith ca += cnzi; 455d50806bdSBarry Smith cj += cnzi; 456d50806bdSBarry Smith } 457716bacf3SKris Buschelman ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 458716bacf3SKris Buschelman ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 459716bacf3SKris Buschelman 460d50806bdSBarry Smith /* Free temp */ 461d50806bdSBarry Smith ierr = PetscFree(temp);CHKERRQ(ierr); 462d50806bdSBarry Smith ierr = PetscLogFlops(flops);CHKERRQ(ierr); 463d50806bdSBarry Smith PetscFunctionReturn(0); 464d50806bdSBarry Smith } 465