1be1d678aSKris Buschelman 22499ec78SHong Zhang /* 32499ec78SHong Zhang Defines matrix-matrix product routines for pairs of MPIAIJ matrices 42499ec78SHong Zhang C = A * B 52499ec78SHong Zhang */ 6c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/aij.h> /*I "petscmat.h" I*/ 7c6db04a5SJed Brown #include <../src/mat/utils/freespace.h> 8c6db04a5SJed Brown #include <../src/mat/impls/aij/mpi/mpiaij.h> 9c6db04a5SJed Brown #include <petscbt.h> 10c6db04a5SJed Brown #include <../src/mat/impls/dense/mpi/mpidense.h> 116291f5a6SHong Zhang /* 126291f5a6SHong Zhang #define DEBUG_MATMATMULT 136291f5a6SHong Zhang */ 141634b032SHong Zhang 152499ec78SHong Zhang #undef __FUNCT__ 162499ec78SHong Zhang #define __FUNCT__ "MatMatMult_MPIAIJ_MPIAIJ" 172499ec78SHong Zhang PetscErrorCode MatMatMult_MPIAIJ_MPIAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill, Mat *C) 182499ec78SHong Zhang { 192499ec78SHong Zhang PetscErrorCode ierr; 202499ec78SHong Zhang 212499ec78SHong Zhang PetscFunctionBegin; 222499ec78SHong Zhang if (scall == MAT_INITIAL_MATRIX){ 23598bc09dSHong Zhang ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr); 24a1a4e84aSHong Zhang ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ(A,B,fill,C);CHKERRQ(ierr); 25598bc09dSHong Zhang ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr); 262499ec78SHong Zhang } 27598bc09dSHong Zhang ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr); 28598bc09dSHong Zhang ierr = (*(*C)->ops->matmultnumeric)(A,B,*C);CHKERRQ(ierr); 29598bc09dSHong Zhang ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr); 302499ec78SHong Zhang PetscFunctionReturn(0); 312499ec78SHong Zhang } 322499ec78SHong Zhang 33f33d1a9aSHong Zhang #undef __FUNCT__ 34776b82aeSLisandro Dalcin #define __FUNCT__ "PetscContainerDestroy_Mat_MatMatMultMPI" 35776b82aeSLisandro Dalcin PetscErrorCode PetscContainerDestroy_Mat_MatMatMultMPI(void *ptr) 36f33d1a9aSHong Zhang { 37f33d1a9aSHong Zhang PetscErrorCode ierr; 389649163dSHong Zhang Mat_MatMatMultMPI *mult=(Mat_MatMatMultMPI*)ptr; 39f33d1a9aSHong Zhang 40f33d1a9aSHong Zhang PetscFunctionBegin; 416bf464f9SBarry Smith ierr = ISDestroy(&mult->isrowa);CHKERRQ(ierr); 426bf464f9SBarry Smith ierr = ISDestroy(&mult->isrowb);CHKERRQ(ierr); 436bf464f9SBarry Smith ierr = ISDestroy(&mult->iscolb);CHKERRQ(ierr); 446bf464f9SBarry Smith ierr = MatDestroy(&mult->C_seq);CHKERRQ(ierr); 456bf464f9SBarry Smith ierr = MatDestroy(&mult->A_loc);CHKERRQ(ierr); 466bf464f9SBarry Smith ierr = MatDestroy(&mult->B_seq);CHKERRQ(ierr); 47f33d1a9aSHong Zhang ierr = PetscFree(mult);CHKERRQ(ierr); 48f33d1a9aSHong Zhang PetscFunctionReturn(0); 49f33d1a9aSHong Zhang } 50f33d1a9aSHong Zhang 51598bc09dSHong Zhang #undef __FUNCT__ 52a1a4e84aSHong Zhang #define __FUNCT__ "MatDestroy_MPIAIJ_MatMatMult" 53a1a4e84aSHong Zhang PetscErrorCode MatDestroy_MPIAIJ_MatMatMult(Mat A) 54598bc09dSHong Zhang { 55598bc09dSHong Zhang PetscErrorCode ierr; 56598bc09dSHong Zhang Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data; 57598bc09dSHong Zhang Mat_PtAPMPI *ptap=a->ptap; 58598bc09dSHong Zhang 59598bc09dSHong Zhang PetscFunctionBegin; 60598bc09dSHong Zhang ierr = PetscFree2(ptap->startsj,ptap->startsj_r);CHKERRQ(ierr); 61598bc09dSHong Zhang ierr = PetscFree(ptap->bufa);CHKERRQ(ierr); 62a1a4e84aSHong Zhang ierr = MatDestroy(&ptap->P_loc);CHKERRQ(ierr); 63a1a4e84aSHong Zhang ierr = MatDestroy(&ptap->P_oth);CHKERRQ(ierr); 64a1a4e84aSHong Zhang ierr = PetscFree(ptap->api);CHKERRQ(ierr); 65a1a4e84aSHong Zhang ierr = PetscFree(ptap->apj);CHKERRQ(ierr); 66598bc09dSHong Zhang ierr = PetscFree(ptap->apa);CHKERRQ(ierr); 67598bc09dSHong Zhang ierr = ptap->destroy(A);CHKERRQ(ierr); 68598bc09dSHong Zhang ierr = PetscFree(ptap);CHKERRQ(ierr); 69598bc09dSHong Zhang PetscFunctionReturn(0); 70598bc09dSHong Zhang } 71598bc09dSHong Zhang 722499ec78SHong Zhang #undef __FUNCT__ 73a1a4e84aSHong Zhang #define __FUNCT__ "MatDestroy_MPIAIJ_MatMatMult_32" 74a1a4e84aSHong Zhang PetscErrorCode MatDestroy_MPIAIJ_MatMatMult_32(Mat A) 752499ec78SHong Zhang { 762499ec78SHong Zhang PetscErrorCode ierr; 77776b82aeSLisandro Dalcin PetscContainer container; 789649163dSHong Zhang Mat_MatMatMultMPI *mult=PETSC_NULL; 792499ec78SHong Zhang 802499ec78SHong Zhang PetscFunctionBegin; 819649163dSHong Zhang ierr = PetscObjectQuery((PetscObject)A,"Mat_MatMatMultMPI",(PetscObject *)&container);CHKERRQ(ierr); 8229825010SHong Zhang if (!container) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_PLIB,"Container does not exist"); 83776b82aeSLisandro Dalcin ierr = PetscContainerGetPointer(container,(void **)&mult);CHKERRQ(ierr); 84dce485f0SBarry Smith A->ops->destroy = mult->destroy; 85bf0cc555SLisandro Dalcin A->ops->duplicate = mult->duplicate; 86bf0cc555SLisandro Dalcin if (A->ops->destroy) { 87992edd73SBarry Smith ierr = (*A->ops->destroy)(A);CHKERRQ(ierr); 88bf0cc555SLisandro Dalcin } 89bf0cc555SLisandro Dalcin ierr = PetscObjectCompose((PetscObject)A,"Mat_MatMatMultMPI",0);CHKERRQ(ierr); 902499ec78SHong Zhang PetscFunctionReturn(0); 912499ec78SHong Zhang } 922499ec78SHong Zhang 932499ec78SHong Zhang #undef __FUNCT__ 94a1a4e84aSHong Zhang #define __FUNCT__ "MatDuplicate_MPIAIJ_MatMatMult_32" 95a1a4e84aSHong Zhang PetscErrorCode MatDuplicate_MPIAIJ_MatMatMult_32(Mat A, MatDuplicateOption op, Mat *M) 96b9c92825SBarry Smith { 97abf897f8SHong Zhang PetscErrorCode ierr; 98abf897f8SHong Zhang Mat_MatMatMultMPI *mult; 99776b82aeSLisandro Dalcin PetscContainer container; 100abf897f8SHong Zhang 101abf897f8SHong Zhang PetscFunctionBegin; 102abf897f8SHong Zhang ierr = PetscObjectQuery((PetscObject)A,"Mat_MatMatMultMPI",(PetscObject *)&container);CHKERRQ(ierr); 10329825010SHong Zhang if (!container) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_PLIB,"Container does not exist"); 104776b82aeSLisandro Dalcin ierr = PetscContainerGetPointer(container,(void **)&mult);CHKERRQ(ierr); 1055c088420SHong Zhang /* Note: the container is not duplicated, because it requires deep copying of 1065c088420SHong Zhang several large data sets (see PetscContainerDestroy_Mat_MatMatMultMPI()). 1075c088420SHong Zhang These data sets are only used for repeated calling of MatMatMultNumeric(). 1085c088420SHong Zhang *M is unlikely being used in this way. Thus we create *M with pure mpiaij format */ 109dce485f0SBarry Smith ierr = (*mult->duplicate)(A,op,M);CHKERRQ(ierr); 110dce485f0SBarry Smith (*M)->ops->destroy = mult->destroy; /* = MatDestroy_MPIAIJ, *M doesn't duplicate A's container! */ 111dce485f0SBarry Smith (*M)->ops->duplicate = mult->duplicate; /* = MatDuplicate_MPIAIJ */ 112abf897f8SHong Zhang PetscFunctionReturn(0); 113abf897f8SHong Zhang } 114abf897f8SHong Zhang 115abf897f8SHong Zhang #undef __FUNCT__ 116a1a4e84aSHong Zhang #define __FUNCT__ "MatDuplicate_MPIAIJ_MatMatMult" 117a1a4e84aSHong Zhang PetscErrorCode MatDuplicate_MPIAIJ_MatMatMult(Mat A, MatDuplicateOption op, Mat *M) 1184ae0847bSHong Zhang { 1194ae0847bSHong Zhang PetscErrorCode ierr; 1204ae0847bSHong Zhang Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data; 1214ae0847bSHong Zhang Mat_PtAPMPI *ptap=a->ptap; 1224ae0847bSHong Zhang 1234ae0847bSHong Zhang PetscFunctionBegin; 1244ae0847bSHong Zhang ierr = (*ptap->duplicate)(A,op,M);CHKERRQ(ierr); 1254ae0847bSHong Zhang (*M)->ops->destroy = ptap->destroy; /* = MatDestroy_MPIAIJ, *M doesn't duplicate A's special structure! */ 1264ae0847bSHong Zhang (*M)->ops->duplicate = ptap->duplicate; /* = MatDuplicate_MPIAIJ */ 1274ae0847bSHong Zhang PetscFunctionReturn(0); 1284ae0847bSHong Zhang } 1294ae0847bSHong Zhang 1304ae0847bSHong Zhang #undef __FUNCT__ 131a1a4e84aSHong Zhang #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIAIJ" 132a1a4e84aSHong Zhang PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C) 133598bc09dSHong Zhang { 134598bc09dSHong Zhang PetscErrorCode ierr; 1354ae0847bSHong Zhang Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data,*c=(Mat_MPIAIJ*)C->data; 136598bc09dSHong Zhang Mat_SeqAIJ *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data; 137598bc09dSHong Zhang Mat_SeqAIJ *cd=(Mat_SeqAIJ*)(c->A)->data,*co=(Mat_SeqAIJ*)(c->B)->data; 138598bc09dSHong Zhang PetscInt *adi=ad->i,*adj,*aoi=ao->i,*aoj; 139598bc09dSHong Zhang PetscScalar *ada,*aoa,*cda=cd->a,*coa=co->a; 140598bc09dSHong Zhang Mat_SeqAIJ *p_loc,*p_oth; 141598bc09dSHong Zhang PetscInt *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pj; 142598bc09dSHong Zhang PetscScalar *pa_loc,*pa_oth,*pa,*apa,valtmp,*ca; 143598bc09dSHong Zhang PetscInt cm=C->rmap->n,anz,pnz; 144598bc09dSHong Zhang Mat_PtAPMPI *ptap=c->ptap; 14565a57a32SSean Farley PetscInt *api,*apj,*apJ,i,j,k,row; 14629825010SHong Zhang PetscInt cstart=C->cmap->rstart; 147598bc09dSHong Zhang PetscInt cdnz,conz,k0,k1; 1487c2e2f26SHong Zhang #if defined(DEBUG_MATMATMULT) 1497c2e2f26SHong Zhang PetscMPIInt rank=a->rank; 1507c2e2f26SHong Zhang #endif 151598bc09dSHong Zhang 152598bc09dSHong Zhang PetscFunctionBegin; 1537c2e2f26SHong Zhang #if defined(DEBUG_MATMATMULT) 1547c2e2f26SHong Zhang if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] call MatMatMultNumeric_MPIAIJ_MPIAIJ()...\n",rank); 1557c2e2f26SHong Zhang #endif 1567c2e2f26SHong Zhang 157a1a4e84aSHong Zhang /* 1) get P_oth = ptap->P_oth and P_loc = ptap->P_loc */ 158598bc09dSHong Zhang /*-----------------------------------------------------*/ 159a1a4e84aSHong Zhang /* update numerical values of P_oth and P_loc */ 160a1a4e84aSHong Zhang ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 161a1a4e84aSHong Zhang ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 1626291f5a6SHong Zhang #if defined(DEBUG_MATMATMULT) 1636291f5a6SHong Zhang if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] got P_oth and P_loc...\n",rank); 1646291f5a6SHong Zhang #endif 165598bc09dSHong Zhang 166598bc09dSHong Zhang /* 2) compute numeric C_loc = A_loc*P = Ad*P_loc + Ao*P_oth */ 167598bc09dSHong Zhang /*----------------------------------------------------------*/ 168598bc09dSHong Zhang /* get data from symbolic products */ 169a1a4e84aSHong Zhang p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 170a1a4e84aSHong Zhang p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 171598bc09dSHong Zhang pi_loc=p_loc->i; pj_loc=p_loc->j; pa_loc=p_loc->a; 172598bc09dSHong Zhang pi_oth=p_oth->i; pj_oth=p_oth->j; pa_oth=p_oth->a; 173598bc09dSHong Zhang 174598bc09dSHong Zhang /* get apa for storing dense row A[i,:]*P */ 175598bc09dSHong Zhang apa = ptap->apa; 176598bc09dSHong Zhang 17729825010SHong Zhang api = ptap->api; 17829825010SHong Zhang apj = ptap->apj; 179598bc09dSHong Zhang for (i=0; i<cm; i++) { 180598bc09dSHong Zhang /* diagonal portion of A */ 181598bc09dSHong Zhang anz = adi[i+1] - adi[i]; 182598bc09dSHong Zhang adj = ad->j + adi[i]; 183598bc09dSHong Zhang ada = ad->a + adi[i]; 184598bc09dSHong Zhang for (j=0; j<anz; j++) { 185598bc09dSHong Zhang row = adj[j]; 186598bc09dSHong Zhang pnz = pi_loc[row+1] - pi_loc[row]; 187598bc09dSHong Zhang pj = pj_loc + pi_loc[row]; 188598bc09dSHong Zhang pa = pa_loc + pi_loc[row]; 189598bc09dSHong Zhang 190598bc09dSHong Zhang /* perform dense axpy */ 191598bc09dSHong Zhang valtmp = ada[j]; 192598bc09dSHong Zhang for (k=0; k<pnz; k++){ 193598bc09dSHong Zhang apa[pj[k]] += valtmp*pa[k]; 194598bc09dSHong Zhang } 195598bc09dSHong Zhang ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr); 196598bc09dSHong Zhang } 197598bc09dSHong Zhang 198598bc09dSHong Zhang /* off-diagonal portion of A */ 199598bc09dSHong Zhang anz = aoi[i+1] - aoi[i]; 200598bc09dSHong Zhang aoj = ao->j + aoi[i]; 201598bc09dSHong Zhang aoa = ao->a + aoi[i]; 202598bc09dSHong Zhang for (j=0; j<anz; j++) { 203598bc09dSHong Zhang row = aoj[j]; 204598bc09dSHong Zhang pnz = pi_oth[row+1] - pi_oth[row]; 205598bc09dSHong Zhang pj = pj_oth + pi_oth[row]; 206598bc09dSHong Zhang pa = pa_oth + pi_oth[row]; 207598bc09dSHong Zhang 208598bc09dSHong Zhang /* perform dense axpy */ 209598bc09dSHong Zhang valtmp = aoa[j]; 210598bc09dSHong Zhang for (k=0; k<pnz; k++){ 211598bc09dSHong Zhang apa[pj[k]] += valtmp*pa[k]; 212598bc09dSHong Zhang } 213598bc09dSHong Zhang ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr); 214598bc09dSHong Zhang } 215598bc09dSHong Zhang 216598bc09dSHong Zhang /* set values in C */ 217598bc09dSHong Zhang apJ = apj + api[i]; 218598bc09dSHong Zhang cdnz = cd->i[i+1] - cd->i[i]; 219598bc09dSHong Zhang conz = co->i[i+1] - co->i[i]; 220598bc09dSHong Zhang 221598bc09dSHong Zhang /* 1st off-diagoanl part of C */ 222598bc09dSHong Zhang ca = coa + co->i[i]; 223598bc09dSHong Zhang k = 0; 224598bc09dSHong Zhang for (k0=0; k0<conz; k0++){ 225598bc09dSHong Zhang if (apJ[k] >= cstart) break; 226598bc09dSHong Zhang ca[k0] = apa[apJ[k]]; 227598bc09dSHong Zhang apa[apJ[k]] = 0.0; 228598bc09dSHong Zhang k++; 229598bc09dSHong Zhang } 230598bc09dSHong Zhang 231598bc09dSHong Zhang /* diagonal part of C */ 232598bc09dSHong Zhang ca = cda + cd->i[i]; 233598bc09dSHong Zhang for (k1=0; k1<cdnz; k1++){ 234598bc09dSHong Zhang ca[k1] = apa[apJ[k]]; 235598bc09dSHong Zhang apa[apJ[k]] = 0.0; 236598bc09dSHong Zhang k++; 237598bc09dSHong Zhang } 238598bc09dSHong Zhang 239598bc09dSHong Zhang /* 2nd off-diagoanl part of C */ 240598bc09dSHong Zhang ca = coa + co->i[i]; 241598bc09dSHong Zhang for (; k0<conz; k0++){ 242598bc09dSHong Zhang ca[k0] = apa[apJ[k]]; 243598bc09dSHong Zhang apa[apJ[k]] = 0.0; 244598bc09dSHong Zhang k++; 245598bc09dSHong Zhang } 246598bc09dSHong Zhang } 247598bc09dSHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 248598bc09dSHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 249598bc09dSHong Zhang PetscFunctionReturn(0); 250598bc09dSHong Zhang } 251598bc09dSHong Zhang 252598bc09dSHong Zhang #undef __FUNCT__ 253a1a4e84aSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ" 254a1a4e84aSHong Zhang PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C) 255598bc09dSHong Zhang { 256598bc09dSHong Zhang PetscErrorCode ierr; 257598bc09dSHong Zhang MPI_Comm comm=((PetscObject)A)->comm; 258bfcd1a73SHong Zhang Mat Cmpi; 259598bc09dSHong Zhang Mat_PtAPMPI *ptap; 260598bc09dSHong Zhang PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 2614ae0847bSHong Zhang Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data,*c; 262bfcd1a73SHong Zhang Mat_SeqAIJ *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_loc,*p_oth; 2634ae0847bSHong Zhang PetscInt *pi_loc,*pj_loc,*pi_oth,*pj_oth,*dnz,*onz; 2644ae0847bSHong Zhang PetscInt *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,rstart=A->rmap->rstart; 265bfcd1a73SHong Zhang PetscInt nlnk,*lnk,i,pnz,row,*api,*apj,*Jptr,apnz,nspacedouble=0,j,nzi; 266bfcd1a73SHong Zhang PetscInt am=A->rmap->n,pN=P->cmap->N,pn=P->cmap->n,pm=P->rmap->n; 267598bc09dSHong Zhang PetscBT lnkbt; 268598bc09dSHong Zhang PetscScalar *apa; 269bfcd1a73SHong Zhang PetscReal afill; 270cf1a0672SHong Zhang PetscBool scalable=PETSC_FALSE; 2717c2e2f26SHong Zhang #if defined(DEBUG_MATMATMULT) 2727c2e2f26SHong Zhang PetscMPIInt rank=a->rank; 2737c2e2f26SHong Zhang #endif 274598bc09dSHong Zhang 275598bc09dSHong Zhang PetscFunctionBegin; 276a1a4e84aSHong Zhang if (A->cmap->rstart != P->rmap->rstart || A->cmap->rend != P->rmap->rend){ 277cf1a0672SHong Zhang SETERRQ4(comm,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, (%D, %D) != (%D,%D)",A->cmap->rstart,A->cmap->rend,P->rmap->rstart,P->rmap->rend); 278a1a4e84aSHong Zhang } 279152983bfSHong Zhang 280152983bfSHong Zhang ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr); 281cf1a0672SHong Zhang 282cf1a0672SHong Zhang ierr = PetscOptionsBool("-matmatmult_32","Use a scalable but slower C=A*B","",scalable,&scalable,PETSC_NULL);CHKERRQ(ierr); 283cf1a0672SHong Zhang if (scalable){ 284cf1a0672SHong Zhang ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable_32(A,P,fill,C);;CHKERRQ(ierr); 285a1a4e84aSHong Zhang PetscFunctionReturn(0); 286a1a4e84aSHong Zhang } 287cf1a0672SHong Zhang ierr = PetscOptionsBool("-matmatmult_scalable","Use a scalable but slower C=A*B","",scalable,&scalable,PETSC_NULL);CHKERRQ(ierr); 288cf1a0672SHong Zhang if (scalable){ 28925023028SHong Zhang ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable(A,P,fill,C);;CHKERRQ(ierr); 2901634b032SHong Zhang PetscFunctionReturn(0); 2911634b032SHong Zhang } 292152983bfSHong Zhang ierr = PetscOptionsEnd();CHKERRQ(ierr); 293a1a4e84aSHong Zhang 2947c2e2f26SHong Zhang #if defined(DEBUG_MATMATMULT) 2957c2e2f26SHong Zhang if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] call MatMatMultSymbolic_MPIAIJ_MPIAIJ()...\n",rank); 2967c2e2f26SHong Zhang #endif 2977c2e2f26SHong Zhang 298598bc09dSHong Zhang /* create struct Mat_PtAPMPI and attached it to C later */ 299598bc09dSHong Zhang ierr = PetscNew(Mat_PtAPMPI,&ptap);CHKERRQ(ierr); 300598bc09dSHong Zhang 301598bc09dSHong Zhang /* malloc apa to store dense row A[i,:]*P */ 302bfcd1a73SHong Zhang ierr = PetscMalloc(pN*sizeof(PetscScalar),&apa);CHKERRQ(ierr); 303bfcd1a73SHong Zhang ierr = PetscMemzero(apa,pN*sizeof(PetscScalar));CHKERRQ(ierr); 304598bc09dSHong Zhang ptap->apa = apa; 3056291f5a6SHong Zhang #if defined(DEBUG_MATMATMULT) 306e9e4536cSHong Zhang if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] Malloc apa pN %D is done...\n",rank,pN); 3076291f5a6SHong Zhang #endif 308598bc09dSHong Zhang 309598bc09dSHong Zhang /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */ 310a1a4e84aSHong Zhang ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 3116291f5a6SHong Zhang #if defined(DEBUG_MATMATMULT) 3126291f5a6SHong Zhang if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] P_oth is done...\n",rank); 3136291f5a6SHong Zhang #endif 314598bc09dSHong Zhang /* get P_loc by taking all local rows of P */ 315a1a4e84aSHong Zhang ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 316598bc09dSHong Zhang 317a1a4e84aSHong Zhang p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 318a1a4e84aSHong Zhang p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 319598bc09dSHong Zhang pi_loc = p_loc->i; pj_loc = p_loc->j; 320598bc09dSHong Zhang pi_oth = p_oth->i; pj_oth = p_oth->j; 321598bc09dSHong Zhang 3226291f5a6SHong Zhang #if defined(DEBUG_MATMATMULT) 323e9e4536cSHong Zhang if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] P_loc is done, Annz %D * P_locnnz %D = %D...\n",rank,ad->i[am]+ao->i[am],p_loc->rmax,(ad->i[am]+ao->i[am])*p_loc->rmax); 3246291f5a6SHong Zhang #endif 3256291f5a6SHong Zhang 326598bc09dSHong Zhang /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */ 327598bc09dSHong Zhang /*-------------------------------------------------------------------*/ 328598bc09dSHong Zhang ierr = PetscMalloc((am+2)*sizeof(PetscInt),&api);CHKERRQ(ierr); 329a1a4e84aSHong Zhang ptap->api = api; 330598bc09dSHong Zhang api[0] = 0; 331598bc09dSHong Zhang 332598bc09dSHong Zhang /* create and initialize a linked list */ 333598bc09dSHong Zhang nlnk = pN+1; 334598bc09dSHong Zhang ierr = PetscLLCreate(pN,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr); 335598bc09dSHong Zhang 336bfcd1a73SHong Zhang /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */ 337bfcd1a73SHong Zhang ierr = PetscFreeSpaceGet((PetscInt)(fill*(adi[am]+aoi[am]+pi_loc[pm])),&free_space);CHKERRQ(ierr); 338598bc09dSHong Zhang current_space = free_space; 339598bc09dSHong Zhang 340598bc09dSHong Zhang ierr = MatPreallocateInitialize(comm,am,pn,dnz,onz);CHKERRQ(ierr); 341598bc09dSHong Zhang for (i=0; i<am; i++) { 342598bc09dSHong Zhang apnz = 0; 343598bc09dSHong Zhang /* diagonal portion of A */ 344598bc09dSHong Zhang nzi = adi[i+1] - adi[i]; 345598bc09dSHong Zhang for (j=0; j<nzi; j++){ 346598bc09dSHong Zhang row = *adj++; 347598bc09dSHong Zhang pnz = pi_loc[row+1] - pi_loc[row]; 348598bc09dSHong Zhang Jptr = pj_loc + pi_loc[row]; 349598bc09dSHong Zhang /* add non-zero cols of P into the sorted linked list lnk */ 350dadf0e6bSHong Zhang ierr = PetscLLAddSorted(pnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr); 351598bc09dSHong Zhang apnz += nlnk; 352598bc09dSHong Zhang } 353598bc09dSHong Zhang /* off-diagonal portion of A */ 354598bc09dSHong Zhang nzi = aoi[i+1] - aoi[i]; 355598bc09dSHong Zhang for (j=0; j<nzi; j++){ 356598bc09dSHong Zhang row = *aoj++; 357598bc09dSHong Zhang pnz = pi_oth[row+1] - pi_oth[row]; 358598bc09dSHong Zhang Jptr = pj_oth + pi_oth[row]; 359dadf0e6bSHong Zhang ierr = PetscLLAddSorted(pnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr); 360598bc09dSHong Zhang apnz += nlnk; 361598bc09dSHong Zhang } 362598bc09dSHong Zhang 363598bc09dSHong Zhang api[i+1] = api[i] + apnz; 364598bc09dSHong Zhang 365598bc09dSHong Zhang /* if free space is not available, double the total space in the list */ 366598bc09dSHong Zhang if (current_space->local_remaining<apnz) { 367598bc09dSHong Zhang ierr = PetscFreeSpaceGet(apnz+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 368598bc09dSHong Zhang nspacedouble++; 369598bc09dSHong Zhang } 370598bc09dSHong Zhang 371598bc09dSHong Zhang /* Copy data into free space, then initialize lnk */ 372598bc09dSHong Zhang ierr = PetscLLClean(pN,pN,apnz,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 373598bc09dSHong Zhang ierr = MatPreallocateSet(i+rstart,apnz,current_space->array,dnz,onz);CHKERRQ(ierr); 374598bc09dSHong Zhang current_space->array += apnz; 375598bc09dSHong Zhang current_space->local_used += apnz; 376598bc09dSHong Zhang current_space->local_remaining -= apnz; 377598bc09dSHong Zhang } 378598bc09dSHong Zhang 379598bc09dSHong Zhang /* Allocate space for apj, initialize apj, and */ 380598bc09dSHong Zhang /* destroy list of free space and other temporary array(s) */ 381a1a4e84aSHong Zhang ierr = PetscMalloc((api[am]+1)*sizeof(PetscInt),&ptap->apj);CHKERRQ(ierr); 382a1a4e84aSHong Zhang apj = ptap->apj; 383a1a4e84aSHong Zhang ierr = PetscFreeSpaceContiguous(&free_space,ptap->apj);CHKERRQ(ierr); 384598bc09dSHong Zhang ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 3856291f5a6SHong Zhang #if defined(DEBUG_MATMATMULT) 3866291f5a6SHong Zhang if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] AP is done...\n",rank); 3876291f5a6SHong Zhang #endif 388598bc09dSHong Zhang 389bfcd1a73SHong Zhang /* create and assemble symbolic parallel matrix Cmpi */ 390598bc09dSHong Zhang /*----------------------------------------------------*/ 391bfcd1a73SHong Zhang ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr); 392bfcd1a73SHong Zhang ierr = MatSetSizes(Cmpi,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 393bfcd1a73SHong Zhang ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr); 394bfcd1a73SHong Zhang ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr); 395598bc09dSHong Zhang ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 396bfcd1a73SHong Zhang ierr = MatSetBlockSize(Cmpi,1);CHKERRQ(ierr); 397598bc09dSHong Zhang for (i=0; i<am; i++){ 398598bc09dSHong Zhang row = i + rstart; 399598bc09dSHong Zhang apnz = api[i+1] - api[i]; 400bfcd1a73SHong Zhang ierr = MatSetValues(Cmpi,1,&row,apnz,apj,apa,INSERT_VALUES);CHKERRQ(ierr); 401598bc09dSHong Zhang apj += apnz; 402598bc09dSHong Zhang } 403bfcd1a73SHong Zhang ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 404bfcd1a73SHong Zhang ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 405598bc09dSHong Zhang 406bfcd1a73SHong Zhang ptap->destroy = Cmpi->ops->destroy; 407bfcd1a73SHong Zhang ptap->duplicate = Cmpi->ops->duplicate; 408bfcd1a73SHong Zhang Cmpi->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ; 409bfcd1a73SHong Zhang Cmpi->ops->destroy = MatDestroy_MPIAIJ_MatMatMult; 410bfcd1a73SHong Zhang Cmpi->ops->duplicate = MatDuplicate_MPIAIJ_MatMatMult; 411598bc09dSHong Zhang 412bfcd1a73SHong Zhang /* attach the supporting struct to Cmpi for reuse */ 413bfcd1a73SHong Zhang c = (Mat_MPIAIJ*)Cmpi->data; 414598bc09dSHong Zhang c->ptap = ptap; 415598bc09dSHong Zhang 416bfcd1a73SHong Zhang *C = Cmpi; 417bfcd1a73SHong Zhang 418bfcd1a73SHong Zhang /* set MatInfo */ 419bfcd1a73SHong Zhang afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]) + 1.e-5; 420bfcd1a73SHong Zhang if (afill < 1.0) afill = 1.0; 421bfcd1a73SHong Zhang Cmpi->info.mallocs = nspacedouble; 422bfcd1a73SHong Zhang Cmpi->info.fill_ratio_given = fill; 423bfcd1a73SHong Zhang Cmpi->info.fill_ratio_needed = afill; 424bfcd1a73SHong Zhang 425bfcd1a73SHong Zhang #if defined(PETSC_USE_INFO) 426bfcd1a73SHong Zhang if (api[am]) { 427bfcd1a73SHong Zhang ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr); 428bfcd1a73SHong Zhang ierr = PetscInfo1(Cmpi,"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr); 429bfcd1a73SHong Zhang } else { 430bfcd1a73SHong Zhang ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr); 431bfcd1a73SHong Zhang } 432bfcd1a73SHong Zhang #endif 433598bc09dSHong Zhang PetscFunctionReturn(0); 434598bc09dSHong Zhang } 435598bc09dSHong Zhang 436a1a4e84aSHong Zhang /* implementation used in PETSc-3.2 */ 437a1a4e84aSHong Zhang /* This routine is called ONLY in the case of reusing previously computed symbolic C */ 438598bc09dSHong Zhang #undef __FUNCT__ 439cf1a0672SHong Zhang #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIAIJ_Scalable_32" 440cf1a0672SHong Zhang PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ_Scalable_32(Mat A,Mat B,Mat C) 441a1a4e84aSHong Zhang { 442a1a4e84aSHong Zhang PetscErrorCode ierr; 443a1a4e84aSHong Zhang Mat *seq; 444a1a4e84aSHong Zhang Mat_MatMatMultMPI *mult; 445a1a4e84aSHong Zhang PetscContainer container; 446a1a4e84aSHong Zhang 447a1a4e84aSHong Zhang PetscFunctionBegin; 44815f4e0f4SHong Zhang ierr = PetscObjectQuery((PetscObject)C,"Mat_MatMatMultMPI",(PetscObject *)&container);CHKERRQ(ierr); 44929825010SHong Zhang if (!container) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_PLIB,"Container does not exist"); 45015f4e0f4SHong Zhang ierr = PetscContainerGetPointer(container,(void **)&mult);CHKERRQ(ierr); 45115f4e0f4SHong Zhang 45215f4e0f4SHong Zhang if (mult->skipNumeric){ /* first numeric product is done during symbolic product */ 45315f4e0f4SHong Zhang mult->skipNumeric = PETSC_FALSE; 45415f4e0f4SHong Zhang PetscFunctionReturn(0); 45515f4e0f4SHong Zhang } 4567c2e2f26SHong Zhang #if defined(DEBUG_MATMATMULT) 4577c2e2f26SHong Zhang PetscMPIInt rank; 458fe615159SHong Zhang MPI_Comm comm = ((PetscObject)C)->comm; 459fe615159SHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 460fe615159SHong Zhang ierr = MPI_Barrier(comm);CHKERRQ(ierr); 461cf1a0672SHong Zhang if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] call MatMatMultNumeric_MPIAIJ_MPIAIJ_Scalable_32()...\n",rank); 4627c2e2f26SHong Zhang #endif 46315f4e0f4SHong Zhang 464a1a4e84aSHong Zhang seq = &mult->B_seq; 465a1a4e84aSHong Zhang ierr = MatGetSubMatrices(B,1,&mult->isrowb,&mult->iscolb,MAT_REUSE_MATRIX,&seq);CHKERRQ(ierr); 466a1a4e84aSHong Zhang mult->B_seq = *seq; 467a1a4e84aSHong Zhang 468a1a4e84aSHong Zhang seq = &mult->A_loc; 469a1a4e84aSHong Zhang ierr = MatGetSubMatrices(A,1,&mult->isrowa,&mult->isrowb,MAT_REUSE_MATRIX,&seq);CHKERRQ(ierr); 470a1a4e84aSHong Zhang mult->A_loc = *seq; 471a1a4e84aSHong Zhang 47225023028SHong Zhang ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(mult->A_loc,mult->B_seq,mult->C_seq);CHKERRQ(ierr); 473a1a4e84aSHong Zhang 474a1a4e84aSHong Zhang ierr = PetscObjectReference((PetscObject)mult->C_seq);CHKERRQ(ierr); 475a1a4e84aSHong Zhang ierr = MatMerge(((PetscObject)A)->comm,mult->C_seq,B->cmap->n,MAT_REUSE_MATRIX,&C);CHKERRQ(ierr); 476a1a4e84aSHong Zhang PetscFunctionReturn(0); 477a1a4e84aSHong Zhang } 478a1a4e84aSHong Zhang 4797c2e2f26SHong Zhang /* numeric product is computed as well */ 480a1a4e84aSHong Zhang #undef __FUNCT__ 481cf1a0672SHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable_32" 482cf1a0672SHong Zhang PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable_32(Mat A,Mat B,PetscReal fill,Mat *C) 4832499ec78SHong Zhang { 4842499ec78SHong Zhang PetscErrorCode ierr; 4852499ec78SHong Zhang Mat_MatMatMultMPI *mult; 486776b82aeSLisandro Dalcin PetscContainer container; 487d5821860SHong Zhang Mat AB,*seq; 488d5821860SHong Zhang Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data; 489d5821860SHong Zhang PetscInt *idx,i,start,ncols,nzA,nzB,*cmap,imark; 4907c2e2f26SHong Zhang #if defined(DEBUG_MATMATMULT) 491e9e4536cSHong Zhang MPI_Comm comm = ((PetscObject)A)->comm; 4927c2e2f26SHong Zhang PetscMPIInt rank=a->rank; 4937c2e2f26SHong Zhang #endif 4942499ec78SHong Zhang 4952499ec78SHong Zhang PetscFunctionBegin; 4967c2e2f26SHong Zhang #if defined(DEBUG_MATMATMULT) 497cf1a0672SHong Zhang if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] call MatMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable_32()...\n",rank); 4987c2e2f26SHong Zhang #endif 499d0f46423SBarry Smith if (A->cmap->rstart != B->rmap->rstart || A->cmap->rend != B->rmap->rend){ 500cf1a0672SHong Zhang SETERRQ4(((PetscObject)A)->comm,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, (%D, %D) != (%D,%D)",A->cmap->rstart,A->cmap->rend,B->rmap->rstart,B->rmap->rend); 5012499ec78SHong Zhang } 502598bc09dSHong Zhang 5032499ec78SHong Zhang ierr = PetscNew(Mat_MatMatMultMPI,&mult);CHKERRQ(ierr); 5042499ec78SHong Zhang 505d5821860SHong Zhang /* get isrowb: nonzero col of A */ 506d5821860SHong Zhang start = A->cmap->rstart; 507d5821860SHong Zhang cmap = a->garray; 508d5821860SHong Zhang nzA = a->A->cmap->n; 509d5821860SHong Zhang nzB = a->B->cmap->n; 510d5821860SHong Zhang ierr = PetscMalloc((nzA+nzB)*sizeof(PetscInt), &idx);CHKERRQ(ierr); 511d5821860SHong Zhang ncols = 0; 512d5821860SHong Zhang for (i=0; i<nzB; i++) { /* row < local row index */ 513d5821860SHong Zhang if (cmap[i] < start) idx[ncols++] = cmap[i]; 514d5821860SHong Zhang else break; 515d5821860SHong Zhang } 516d5821860SHong Zhang imark = i; 517d5821860SHong Zhang for (i=0; i<nzA; i++) idx[ncols++] = start + i; /* local rows */ 518d5821860SHong Zhang for (i=imark; i<nzB; i++) idx[ncols++] = cmap[i]; /* row > local row index */ 519d5821860SHong Zhang ierr = ISCreateGeneral(PETSC_COMM_SELF,ncols,idx,PETSC_OWN_POINTER,&mult->isrowb);CHKERRQ(ierr); 520d5821860SHong Zhang ierr = ISCreateStride(PETSC_COMM_SELF,B->cmap->N,0,1,&mult->iscolb);CHKERRQ(ierr); 521d5821860SHong Zhang 522d5821860SHong Zhang /* get isrowa: all local rows of A */ 523d5821860SHong Zhang ierr = ISCreateStride(PETSC_COMM_SELF,A->rmap->n,A->rmap->rstart,1,&mult->isrowa);CHKERRQ(ierr); 524d5821860SHong Zhang 525d5821860SHong Zhang /* Below should go to MatMatMultNumeric_MPIAIJ_MPIAIJ() - How to generate C there? */ 5262499ec78SHong Zhang /* create a seq matrix B_seq = submatrix of B by taking rows of B that equal to nonzero col of A */ 527d5821860SHong Zhang ierr = MatGetSubMatrices(B,1,&mult->isrowb,&mult->iscolb,MAT_INITIAL_MATRIX,&seq);CHKERRQ(ierr); 528d5821860SHong Zhang mult->B_seq = *seq; 529d5821860SHong Zhang ierr = PetscFree(seq);CHKERRQ(ierr); 5306291f5a6SHong Zhang #if defined(DEBUG_MATMATMULT) 5316291f5a6SHong Zhang if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] B_seq is done...\n",rank); 5326291f5a6SHong Zhang #endif 5332499ec78SHong Zhang 5342499ec78SHong Zhang /* create a seq matrix A_seq = submatrix of A by taking all local rows of A */ 535d5821860SHong Zhang ierr = MatGetSubMatrices(A,1,&mult->isrowa,&mult->isrowb,MAT_INITIAL_MATRIX,&seq);CHKERRQ(ierr); 536d5821860SHong Zhang mult->A_loc = *seq; 537d5821860SHong Zhang ierr = PetscFree(seq);CHKERRQ(ierr); 5386291f5a6SHong Zhang #if defined(DEBUG_MATMATMULT) 539e9e4536cSHong Zhang if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] A_loc is done...\n",rank); 5406291f5a6SHong Zhang #endif 5412499ec78SHong Zhang 5422499ec78SHong Zhang /* compute C_seq = A_seq * B_seq */ 54325023028SHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(mult->A_loc,mult->B_seq,fill,&mult->C_seq);CHKERRQ(ierr); 5446291f5a6SHong Zhang #if defined(DEBUG_MATMATMULT) 545e9e4536cSHong Zhang ierr = MPI_Barrier(comm);CHKERRQ(ierr); 546e9e4536cSHong Zhang if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] C_seq Symbolic is done...\n",rank); 547e9e4536cSHong Zhang #endif 54825023028SHong Zhang ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(mult->A_loc,mult->B_seq,mult->C_seq);CHKERRQ(ierr); 549e9e4536cSHong Zhang #if defined(DEBUG_MATMATMULT) 550e9e4536cSHong Zhang if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] C_seq Numeric is done...\n",rank); 5516291f5a6SHong Zhang #endif 5522499ec78SHong Zhang 5532499ec78SHong Zhang /* create mpi matrix C by concatinating C_seq */ 5542499ec78SHong Zhang ierr = PetscObjectReference((PetscObject)mult->C_seq);CHKERRQ(ierr); /* prevent C_seq being destroyed by MatMerge() */ 555*9b8102ccSHong Zhang ierr = MatMergeSymbolic(((PetscObject)A)->comm,mult->C_seq,B->cmap->n,&AB);CHKERRQ(ierr); 556*9b8102ccSHong Zhang ierr = MatMergeNumeric(((PetscObject)A)->comm,mult->C_seq,B->cmap->n,AB);CHKERRQ(ierr); 5576291f5a6SHong Zhang #if defined(DEBUG_MATMATMULT) 5586291f5a6SHong Zhang if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] Merge is done...\n",rank); 5596291f5a6SHong Zhang #endif 5602499ec78SHong Zhang 5612499ec78SHong Zhang /* attach the supporting struct to C for reuse of symbolic C */ 562776b82aeSLisandro Dalcin ierr = PetscContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr); 563776b82aeSLisandro Dalcin ierr = PetscContainerSetPointer(container,mult);CHKERRQ(ierr); 564776b82aeSLisandro Dalcin ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_Mat_MatMatMultMPI);CHKERRQ(ierr); 565d5821860SHong Zhang ierr = PetscObjectCompose((PetscObject)AB,"Mat_MatMatMultMPI",(PetscObject)container);CHKERRQ(ierr); 566b160f555SLisandro Dalcin ierr = PetscContainerDestroy(&container);CHKERRQ(ierr); 5673f7ae0dbSHong Zhang mult->skipNumeric = PETSC_TRUE; /* a numeric product is done here */ 568d5821860SHong Zhang mult->destroy = AB->ops->destroy; 569d5821860SHong Zhang mult->duplicate = AB->ops->duplicate; 570cf1a0672SHong Zhang AB->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ_Scalable_32; 571a1a4e84aSHong Zhang AB->ops->destroy = MatDestroy_MPIAIJ_MatMatMult_32; 572a1a4e84aSHong Zhang AB->ops->duplicate = MatDuplicate_MPIAIJ_MatMatMult_32; 5738cdbd757SHong Zhang AB->ops->matmult = MatMatMult_MPIAIJ_MPIAIJ; 574d5821860SHong Zhang *C = AB; 5752499ec78SHong Zhang PetscFunctionReturn(0); 5762499ec78SHong Zhang } 5772499ec78SHong Zhang 5789123193aSHong Zhang #undef __FUNCT__ 5799123193aSHong Zhang #define __FUNCT__ "MatMatMult_MPIAIJ_MPIDense" 5809123193aSHong Zhang PetscErrorCode MatMatMult_MPIAIJ_MPIDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 5819123193aSHong Zhang { 5829123193aSHong Zhang PetscErrorCode ierr; 5839123193aSHong Zhang 5849123193aSHong Zhang PetscFunctionBegin; 5859123193aSHong Zhang if (scall == MAT_INITIAL_MATRIX){ 5869123193aSHong Zhang ierr = MatMatMultSymbolic_MPIAIJ_MPIDense(A,B,fill,C);CHKERRQ(ierr); 5879123193aSHong Zhang } 5889123193aSHong Zhang ierr = MatMatMultNumeric_MPIAIJ_MPIDense(A,B,*C);CHKERRQ(ierr); 5899123193aSHong Zhang PetscFunctionReturn(0); 5909123193aSHong Zhang } 5919123193aSHong Zhang 5924324174eSBarry Smith typedef struct { 5934324174eSBarry Smith Mat workB; 5944324174eSBarry Smith PetscScalar *rvalues,*svalues; 5954324174eSBarry Smith MPI_Request *rwaits,*swaits; 5964324174eSBarry Smith } MPIAIJ_MPIDense; 5974324174eSBarry Smith 5987af9e4fdSHong Zhang #undef __FUNCT__ 5997af9e4fdSHong Zhang #define __FUNCT__ "MPIAIJ_MPIDenseDestroy" 6004324174eSBarry Smith PetscErrorCode MPIAIJ_MPIDenseDestroy(void *ctx) 6014324174eSBarry Smith { 6024324174eSBarry Smith MPIAIJ_MPIDense *contents = (MPIAIJ_MPIDense*) ctx; 6034324174eSBarry Smith PetscErrorCode ierr; 6044324174eSBarry Smith 6054324174eSBarry Smith PetscFunctionBegin; 6066bf464f9SBarry Smith ierr = MatDestroy(&contents->workB);CHKERRQ(ierr); 6074324174eSBarry Smith ierr = PetscFree4(contents->rvalues,contents->svalues,contents->rwaits,contents->swaits);CHKERRQ(ierr); 6084324174eSBarry Smith ierr = PetscFree(contents);CHKERRQ(ierr); 6094324174eSBarry Smith PetscFunctionReturn(0); 6104324174eSBarry Smith } 6114324174eSBarry Smith 6129123193aSHong Zhang #undef __FUNCT__ 6139123193aSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIDense" 6149123193aSHong Zhang PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIDense(Mat A,Mat B,PetscReal fill,Mat *C) 6159123193aSHong Zhang { 6169123193aSHong Zhang PetscErrorCode ierr; 617f32d5d43SBarry Smith Mat_MPIAIJ *aij = (Mat_MPIAIJ*) A->data; 618d0f46423SBarry Smith PetscInt nz = aij->B->cmap->n; 619bf0cc555SLisandro Dalcin PetscContainer container; 6204324174eSBarry Smith MPIAIJ_MPIDense *contents; 6214324174eSBarry Smith VecScatter ctx = aij->Mvctx; 6224324174eSBarry Smith VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata; 6234324174eSBarry Smith VecScatter_MPI_General *to = ( VecScatter_MPI_General*) ctx->todata; 624d0f46423SBarry Smith PetscInt m=A->rmap->n,n=B->cmap->n; 6259123193aSHong Zhang 6269123193aSHong Zhang PetscFunctionBegin; 627cb2480eaSBarry Smith ierr = MatCreate(((PetscObject)B)->comm,C);CHKERRQ(ierr); 628d0f46423SBarry Smith ierr = MatSetSizes(*C,m,n,A->rmap->N,B->cmap->N);CHKERRQ(ierr); 629cb2480eaSBarry Smith ierr = MatSetType(*C,MATMPIDENSE);CHKERRQ(ierr); 630cb2480eaSBarry Smith ierr = MatAssemblyBegin(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 631cb2480eaSBarry Smith ierr = MatAssemblyEnd(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 6328cdbd757SHong Zhang (*C)->ops->matmult = MatMatMult_MPIAIJ_MPIDense; 633f32d5d43SBarry Smith 6344324174eSBarry Smith ierr = PetscNew(MPIAIJ_MPIDense,&contents);CHKERRQ(ierr); 635f32d5d43SBarry Smith /* Create work matrix used to store off processor rows of B needed for local product */ 636d0f46423SBarry Smith ierr = MatCreateSeqDense(PETSC_COMM_SELF,nz,B->cmap->N,PETSC_NULL,&contents->workB);CHKERRQ(ierr); 6374324174eSBarry Smith /* Create work arrays needed */ 638d0f46423SBarry Smith ierr = PetscMalloc4(B->cmap->N*from->starts[from->n],PetscScalar,&contents->rvalues, 639d0f46423SBarry Smith B->cmap->N*to->starts[to->n],PetscScalar,&contents->svalues, 6404324174eSBarry Smith from->n,MPI_Request,&contents->rwaits, 6414324174eSBarry Smith to->n,MPI_Request,&contents->swaits);CHKERRQ(ierr); 6424324174eSBarry Smith 643bf0cc555SLisandro Dalcin ierr = PetscContainerCreate(((PetscObject)A)->comm,&container);CHKERRQ(ierr); 644bf0cc555SLisandro Dalcin ierr = PetscContainerSetPointer(container,contents);CHKERRQ(ierr); 645bf0cc555SLisandro Dalcin ierr = PetscContainerSetUserDestroy(container,MPIAIJ_MPIDenseDestroy);CHKERRQ(ierr); 646bf0cc555SLisandro Dalcin ierr = PetscObjectCompose((PetscObject)(*C),"workB",(PetscObject)container);CHKERRQ(ierr); 647bf0cc555SLisandro Dalcin ierr = PetscContainerDestroy(&container);CHKERRQ(ierr); 6489123193aSHong Zhang PetscFunctionReturn(0); 6499123193aSHong Zhang } 6509123193aSHong Zhang 6517af9e4fdSHong Zhang #undef __FUNCT__ 6527af9e4fdSHong Zhang #define __FUNCT__ "MatMPIDenseScatter" 653f32d5d43SBarry Smith /* 6542636ff26SBarry Smith Performs an efficient scatter on the rows of B needed by this process; this is 6552636ff26SBarry Smith a modification of the VecScatterBegin_() routines. 656f32d5d43SBarry Smith */ 6574324174eSBarry Smith PetscErrorCode MatMPIDenseScatter(Mat A,Mat B,Mat C,Mat *outworkB) 658f32d5d43SBarry Smith { 659f32d5d43SBarry Smith Mat_MPIAIJ *aij = (Mat_MPIAIJ*)A->data; 660f32d5d43SBarry Smith PetscErrorCode ierr; 661f32d5d43SBarry Smith PetscScalar *b,*w,*svalues,*rvalues; 662f32d5d43SBarry Smith VecScatter ctx = aij->Mvctx; 663f32d5d43SBarry Smith VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata; 664f32d5d43SBarry Smith VecScatter_MPI_General *to = ( VecScatter_MPI_General*) ctx->todata; 665f32d5d43SBarry Smith PetscInt i,j,k; 666aa5bb8c0SSatish Balay PetscInt *sindices,*sstarts,*rindices,*rstarts; 667aa5bb8c0SSatish Balay PetscMPIInt *sprocs,*rprocs,nrecvs; 668f32d5d43SBarry Smith MPI_Request *swaits,*rwaits; 6697adad957SLisandro Dalcin MPI_Comm comm = ((PetscObject)A)->comm; 670d0f46423SBarry Smith PetscMPIInt tag = ((PetscObject)ctx)->tag,ncols = B->cmap->N, nrows = aij->B->cmap->n,imdex,nrowsB = B->rmap->n; 671f32d5d43SBarry Smith MPI_Status status; 6724324174eSBarry Smith MPIAIJ_MPIDense *contents; 673bf0cc555SLisandro Dalcin PetscContainer container; 6744324174eSBarry Smith Mat workB; 675f32d5d43SBarry Smith 676f32d5d43SBarry Smith PetscFunctionBegin; 677bf0cc555SLisandro Dalcin ierr = PetscObjectQuery((PetscObject)C,"workB",(PetscObject*)&container);CHKERRQ(ierr); 67829825010SHong Zhang if (!container) SETERRQ(comm,PETSC_ERR_PLIB,"Container does not exist"); 679bf0cc555SLisandro Dalcin ierr = PetscContainerGetPointer(container,(void**)&contents);CHKERRQ(ierr); 6804324174eSBarry Smith 6814324174eSBarry Smith workB = *outworkB = contents->workB; 682cf1a0672SHong Zhang if (nrows != workB->rmap->n) SETERRQ2(comm,PETSC_ERR_PLIB,"Number of rows of workB %D not equal to columns of aij->B %D",nrows,workB->cmap->n); 683f32d5d43SBarry Smith sindices = to->indices; 684f32d5d43SBarry Smith sstarts = to->starts; 685f32d5d43SBarry Smith sprocs = to->procs; 6864324174eSBarry Smith swaits = contents->swaits; 6874324174eSBarry Smith svalues = contents->svalues; 688f32d5d43SBarry Smith 689f32d5d43SBarry Smith rindices = from->indices; 690f32d5d43SBarry Smith rstarts = from->starts; 691f32d5d43SBarry Smith rprocs = from->procs; 6924324174eSBarry Smith rwaits = contents->rwaits; 6934324174eSBarry Smith rvalues = contents->rvalues; 694f32d5d43SBarry Smith 695f32d5d43SBarry Smith ierr = MatGetArray(B,&b);CHKERRQ(ierr); 696f32d5d43SBarry Smith ierr = MatGetArray(workB,&w);CHKERRQ(ierr); 697f32d5d43SBarry Smith 698f32d5d43SBarry Smith for (i=0; i<from->n; i++) { 699f32d5d43SBarry Smith ierr = MPI_Irecv(rvalues+ncols*rstarts[i],ncols*(rstarts[i+1]-rstarts[i]),MPIU_SCALAR,rprocs[i],tag,comm,rwaits+i);CHKERRQ(ierr); 700f32d5d43SBarry Smith } 701f32d5d43SBarry Smith 702f32d5d43SBarry Smith for (i=0; i<to->n; i++) { 703f32d5d43SBarry Smith /* pack a message at a time */ 704f32d5d43SBarry Smith CHKMEMQ; 705f32d5d43SBarry Smith for (j=0; j<sstarts[i+1]-sstarts[i]; j++){ 706f32d5d43SBarry Smith for (k=0; k<ncols; k++) { 707f32d5d43SBarry Smith svalues[ncols*(sstarts[i] + j) + k] = b[sindices[sstarts[i]+j] + nrowsB*k]; 7082636ff26SBarry Smith } 7092636ff26SBarry Smith } 710f32d5d43SBarry Smith CHKMEMQ; 711f32d5d43SBarry Smith ierr = MPI_Isend(svalues+ncols*sstarts[i],ncols*(sstarts[i+1]-sstarts[i]),MPIU_SCALAR,sprocs[i],tag,comm,swaits+i);CHKERRQ(ierr); 712f32d5d43SBarry Smith } 7132636ff26SBarry Smith 714f32d5d43SBarry Smith nrecvs = from->n; 715f32d5d43SBarry Smith while (nrecvs) { 716f32d5d43SBarry Smith ierr = MPI_Waitany(from->n,rwaits,&imdex,&status);CHKERRQ(ierr); 717f32d5d43SBarry Smith nrecvs--; 718f32d5d43SBarry Smith /* unpack a message at a time */ 719f32d5d43SBarry Smith CHKMEMQ; 720f32d5d43SBarry Smith for (j=0; j<rstarts[imdex+1]-rstarts[imdex]; j++){ 721f32d5d43SBarry Smith for (k=0; k<ncols; k++) { 722f32d5d43SBarry Smith w[rindices[rstarts[imdex]+j] + nrows*k] = rvalues[ncols*(rstarts[imdex] + j) + k]; 7232636ff26SBarry Smith } 7242636ff26SBarry Smith } 725f32d5d43SBarry Smith CHKMEMQ; 726f32d5d43SBarry Smith } 727cb9801acSJed Brown if (to->n) {ierr = MPI_Waitall(to->n,swaits,to->sstatus);CHKERRQ(ierr);} 728f32d5d43SBarry Smith 729f32d5d43SBarry Smith ierr = MatRestoreArray(B,&b);CHKERRQ(ierr); 730f32d5d43SBarry Smith ierr = MatRestoreArray(workB,&w);CHKERRQ(ierr); 731f32d5d43SBarry Smith ierr = MatAssemblyBegin(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 732f32d5d43SBarry Smith ierr = MatAssemblyEnd(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 733f32d5d43SBarry Smith PetscFunctionReturn(0); 734f32d5d43SBarry Smith } 735f32d5d43SBarry Smith extern PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat,Mat,Mat); 736f32d5d43SBarry Smith 7379123193aSHong Zhang #undef __FUNCT__ 7389123193aSHong Zhang #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIDense" 7399123193aSHong Zhang PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIDense(Mat A,Mat B,Mat C) 7409123193aSHong Zhang { 7419123193aSHong Zhang PetscErrorCode ierr; 742f32d5d43SBarry Smith Mat_MPIAIJ *aij = (Mat_MPIAIJ*)A->data; 743f32d5d43SBarry Smith Mat_MPIDense *bdense = (Mat_MPIDense*)B->data; 744f32d5d43SBarry Smith Mat_MPIDense *cdense = (Mat_MPIDense*)C->data; 745f32d5d43SBarry Smith Mat workB; 7469123193aSHong Zhang 7479123193aSHong Zhang PetscFunctionBegin; 7489123193aSHong Zhang 749f32d5d43SBarry Smith /* diagonal block of A times all local rows of B*/ 750f32d5d43SBarry Smith ierr = MatMatMultNumeric_SeqAIJ_SeqDense(aij->A,bdense->A,cdense->A);CHKERRQ(ierr); 751f32d5d43SBarry Smith 752f32d5d43SBarry Smith /* get off processor parts of B needed to complete the product */ 7534324174eSBarry Smith ierr = MatMPIDenseScatter(A,B,C,&workB);CHKERRQ(ierr); 754f32d5d43SBarry Smith 755f32d5d43SBarry Smith /* off-diagonal block of A times nonlocal rows of B */ 756f32d5d43SBarry Smith ierr = MatMatMultNumericAdd_SeqAIJ_SeqDense(aij->B,workB,cdense->A);CHKERRQ(ierr); 7579123193aSHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 7589123193aSHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 7599123193aSHong Zhang PetscFunctionReturn(0); 7609123193aSHong Zhang } 761cf1a0672SHong Zhang 7621634b032SHong Zhang #undef __FUNCT__ 763cf1a0672SHong Zhang #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIAIJ_Scalable" 764cf1a0672SHong Zhang PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ_Scalable(Mat A,Mat P,Mat C) 7651634b032SHong Zhang { 7661634b032SHong Zhang PetscErrorCode ierr; 767cf1a0672SHong Zhang Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data,*c=(Mat_MPIAIJ*)C->data; 768cf1a0672SHong Zhang Mat_SeqAIJ *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data; 769cf1a0672SHong Zhang Mat_SeqAIJ *cd=(Mat_SeqAIJ*)(c->A)->data,*co=(Mat_SeqAIJ*)(c->B)->data; 770cf1a0672SHong Zhang PetscInt *adi=ad->i,*adj,*aoi=ao->i,*aoj; 771cf1a0672SHong Zhang PetscScalar *ada,*aoa,*cda=cd->a,*coa=co->a; 772cf1a0672SHong Zhang Mat_SeqAIJ *p_loc,*p_oth; 773cf1a0672SHong Zhang PetscInt *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pj; 77429825010SHong Zhang PetscScalar *pa_loc,*pa_oth,*pa,valtmp,*ca; 77529825010SHong Zhang PetscInt cm=C->rmap->n,anz,pnz; 776cf1a0672SHong Zhang Mat_PtAPMPI *ptap=c->ptap; 77729825010SHong Zhang PetscScalar *apa_sparse=ptap->apa; 778cf1a0672SHong Zhang PetscInt *api,*apj,*apJ,i,j,k,row; 77929825010SHong Zhang PetscInt cstart=C->cmap->rstart; 78029825010SHong Zhang PetscInt cdnz,conz,k0,k1,nextp; 7811634b032SHong Zhang #if defined(DEBUG_MATMATMULT) 782cf1a0672SHong Zhang PetscMPIInt rank=a->rank; 7831634b032SHong Zhang #endif 7841634b032SHong Zhang 7851634b032SHong Zhang PetscFunctionBegin; 7861634b032SHong Zhang #if defined(DEBUG_MATMATMULT) 787cf1a0672SHong Zhang if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] call MatMatMultNumeric_MPIAIJ_MPIAIJ_Scalable()...\n",rank); 7881634b032SHong Zhang #endif 789cf1a0672SHong Zhang 790cf1a0672SHong Zhang /* 1) get P_oth = ptap->P_oth and P_loc = ptap->P_loc */ 791cf1a0672SHong Zhang /*-----------------------------------------------------*/ 792cf1a0672SHong Zhang /* update numerical values of P_oth and P_loc */ 793cf1a0672SHong Zhang ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 794cf1a0672SHong Zhang ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 795cf1a0672SHong Zhang #if defined(DEBUG_MATMATMULT) 796cf1a0672SHong Zhang if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] got P_oth and P_loc...\n",rank); 797cf1a0672SHong Zhang #endif 798cf1a0672SHong Zhang 799cf1a0672SHong Zhang /* 2) compute numeric C_loc = A_loc*P = Ad*P_loc + Ao*P_oth */ 800cf1a0672SHong Zhang /*----------------------------------------------------------*/ 801cf1a0672SHong Zhang /* get data from symbolic products */ 802cf1a0672SHong Zhang p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 803cf1a0672SHong Zhang p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 804cf1a0672SHong Zhang pi_loc=p_loc->i; pj_loc=p_loc->j; pa_loc=p_loc->a; 805cf1a0672SHong Zhang pi_oth=p_oth->i; pj_oth=p_oth->j; pa_oth=p_oth->a; 806cf1a0672SHong Zhang 80729825010SHong Zhang api = ptap->api; 80829825010SHong Zhang apj = ptap->apj; 809cf1a0672SHong Zhang for (i=0; i<cm; i++) { 81029825010SHong Zhang apJ = apj + api[i]; 81129825010SHong Zhang 812cf1a0672SHong Zhang /* diagonal portion of A */ 813cf1a0672SHong Zhang anz = adi[i+1] - adi[i]; 814cf1a0672SHong Zhang adj = ad->j + adi[i]; 815cf1a0672SHong Zhang ada = ad->a + adi[i]; 816cf1a0672SHong Zhang for (j=0; j<anz; j++) { 817cf1a0672SHong Zhang row = adj[j]; 818cf1a0672SHong Zhang pnz = pi_loc[row+1] - pi_loc[row]; 819cf1a0672SHong Zhang pj = pj_loc + pi_loc[row]; 820cf1a0672SHong Zhang pa = pa_loc + pi_loc[row]; 82129825010SHong Zhang /* perform sparse axpy */ 822cf1a0672SHong Zhang valtmp = ada[j]; 82329825010SHong Zhang nextp = 0; 82429825010SHong Zhang for (k=0; nextp<pnz; k++) { 82529825010SHong Zhang if (apJ[k] == pj[nextp]) { /* column of AP == column of P */ 82629825010SHong Zhang apa_sparse[k] += valtmp*pa[nextp++]; 82729825010SHong Zhang } 8281634b032SHong Zhang } 829cf1a0672SHong Zhang ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr); 830cf1a0672SHong Zhang } 8311634b032SHong Zhang 832cf1a0672SHong Zhang /* off-diagonal portion of A */ 833cf1a0672SHong Zhang anz = aoi[i+1] - aoi[i]; 834cf1a0672SHong Zhang aoj = ao->j + aoi[i]; 835cf1a0672SHong Zhang aoa = ao->a + aoi[i]; 836cf1a0672SHong Zhang for (j=0; j<anz; j++) { 837cf1a0672SHong Zhang row = aoj[j]; 838cf1a0672SHong Zhang pnz = pi_oth[row+1] - pi_oth[row]; 839cf1a0672SHong Zhang pj = pj_oth + pi_oth[row]; 840cf1a0672SHong Zhang pa = pa_oth + pi_oth[row]; 84129825010SHong Zhang /* perform sparse axpy */ 842cf1a0672SHong Zhang valtmp = aoa[j]; 84329825010SHong Zhang nextp = 0; 84429825010SHong Zhang for (k=0; nextp<pnz; k++) { 84529825010SHong Zhang if (apJ[k] == pj[nextp]) { /* column of AP == column of P */ 84629825010SHong Zhang apa_sparse[k] += valtmp*pa[nextp++]; 84729825010SHong Zhang } 848cf1a0672SHong Zhang } 849cf1a0672SHong Zhang ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr); 850cf1a0672SHong Zhang } 8511634b032SHong Zhang 852cf1a0672SHong Zhang /* set values in C */ 853cf1a0672SHong Zhang cdnz = cd->i[i+1] - cd->i[i]; 854cf1a0672SHong Zhang conz = co->i[i+1] - co->i[i]; 8551634b032SHong Zhang 856cf1a0672SHong Zhang /* 1st off-diagoanl part of C */ 857cf1a0672SHong Zhang ca = coa + co->i[i]; 858cf1a0672SHong Zhang k = 0; 859cf1a0672SHong Zhang for (k0=0; k0<conz; k0++){ 860cf1a0672SHong Zhang if (apJ[k] >= cstart) break; 86129825010SHong Zhang ca[k0] = apa_sparse[k]; 86229825010SHong Zhang apa_sparse[k] = 0.0; 863cf1a0672SHong Zhang k++; 864cf1a0672SHong Zhang } 8651634b032SHong Zhang 866cf1a0672SHong Zhang /* diagonal part of C */ 867cf1a0672SHong Zhang ca = cda + cd->i[i]; 868cf1a0672SHong Zhang for (k1=0; k1<cdnz; k1++){ 86929825010SHong Zhang ca[k1] = apa_sparse[k]; 87029825010SHong Zhang apa_sparse[k] = 0.0; 871cf1a0672SHong Zhang k++; 872cf1a0672SHong Zhang } 873cf1a0672SHong Zhang 874cf1a0672SHong Zhang /* 2nd off-diagoanl part of C */ 875cf1a0672SHong Zhang ca = coa + co->i[i]; 876cf1a0672SHong Zhang for (; k0<conz; k0++){ 87729825010SHong Zhang ca[k0] = apa_sparse[k]; 87829825010SHong Zhang apa_sparse[k] = 0.0; 879cf1a0672SHong Zhang k++; 880cf1a0672SHong Zhang } 881cf1a0672SHong Zhang } 882cf1a0672SHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 883cf1a0672SHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 884cf1a0672SHong Zhang PetscFunctionReturn(0); 885cf1a0672SHong Zhang } 886cf1a0672SHong Zhang 887cf1a0672SHong Zhang /* same as MatMatMultSymbolic_MPIAIJ_MPIAIJ(), except using LLCondensed to avoid O(BN) memory requirement */ 888cf1a0672SHong Zhang #undef __FUNCT__ 889cf1a0672SHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable" 890cf1a0672SHong Zhang PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable(Mat A,Mat P,PetscReal fill,Mat *C) 891cf1a0672SHong Zhang { 892cf1a0672SHong Zhang PetscErrorCode ierr; 893cf1a0672SHong Zhang MPI_Comm comm=((PetscObject)A)->comm; 894cf1a0672SHong Zhang Mat Cmpi; 895cf1a0672SHong Zhang Mat_PtAPMPI *ptap; 896cf1a0672SHong Zhang PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 897cf1a0672SHong Zhang Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data,*c; 898cf1a0672SHong Zhang Mat_SeqAIJ *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_loc,*p_oth; 899cf1a0672SHong Zhang PetscInt *pi_loc,*pj_loc,*pi_oth,*pj_oth,*dnz,*onz; 900cf1a0672SHong Zhang PetscInt *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,rstart=A->rmap->rstart; 901cf1a0672SHong Zhang PetscInt i,pnz,row,*api,*apj,*Jptr,apnz,nspacedouble=0,j,nzi,*nlnk,*lnk,apnz_max=0; 902cf1a0672SHong Zhang PetscInt am=A->rmap->n,pN=P->cmap->N,pn=P->cmap->n,pm=P->rmap->n; 903cf1a0672SHong Zhang PetscInt nlnk_max,armax,prmax; 904cf1a0672SHong Zhang PetscBT lnkbt; 905cf1a0672SHong Zhang PetscReal afill; 906cf1a0672SHong Zhang PetscScalar *apa; 907cf1a0672SHong Zhang #if defined(DEBUG_MATMATMULT) 908cf1a0672SHong Zhang PetscMPIInt rank=a->rank; 909cf1a0672SHong Zhang #endif 910cf1a0672SHong Zhang 911cf1a0672SHong Zhang PetscFunctionBegin; 912cf1a0672SHong Zhang #if defined(DEBUG_MATMATMULT) 913cf1a0672SHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 914cf1a0672SHong Zhang if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] call MatMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable()...\n",rank); 915cf1a0672SHong Zhang #endif 916cf1a0672SHong Zhang 917cf1a0672SHong Zhang /* create struct Mat_PtAPMPI and attached it to C later */ 918cf1a0672SHong Zhang ierr = PetscNew(Mat_PtAPMPI,&ptap);CHKERRQ(ierr); 919cf1a0672SHong Zhang 920cf1a0672SHong Zhang /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */ 921cf1a0672SHong Zhang ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 922cf1a0672SHong Zhang #if defined(DEBUG_MATMATMULT) 923cf1a0672SHong Zhang if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] P_oth is done...\n",rank); 924cf1a0672SHong Zhang #endif 925cf1a0672SHong Zhang /* get P_loc by taking all local rows of P */ 926cf1a0672SHong Zhang ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 927cf1a0672SHong Zhang 928cf1a0672SHong Zhang p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 929cf1a0672SHong Zhang p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 930cf1a0672SHong Zhang pi_loc = p_loc->i; pj_loc = p_loc->j; 931cf1a0672SHong Zhang pi_oth = p_oth->i; pj_oth = p_oth->j; 932cf1a0672SHong Zhang 933cf1a0672SHong Zhang #if defined(DEBUG_MATMATMULT) 934cf1a0672SHong Zhang if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] P_loc is done, Annz %D * P_locnnz %D = %D...\n",rank,ad->i[am]+ao->i[am],p_loc->rmax,(ad->i[am]+ao->i[am])*p_loc->rmax); 935cf1a0672SHong Zhang #endif 936cf1a0672SHong Zhang 937cf1a0672SHong Zhang /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */ 938cf1a0672SHong Zhang /*-------------------------------------------------------------------*/ 939cf1a0672SHong Zhang ierr = PetscMalloc((am+2)*sizeof(PetscInt),&api);CHKERRQ(ierr); 940cf1a0672SHong Zhang ptap->api = api; 941cf1a0672SHong Zhang api[0] = 0; 942cf1a0672SHong Zhang 943cf1a0672SHong Zhang /* create and initialize a linked list */ 944cf1a0672SHong Zhang armax = ad->rmax+ao->rmax; 945cf1a0672SHong Zhang prmax = PetscMax(p_loc->rmax,p_oth->rmax); 946cf1a0672SHong Zhang nlnk_max = armax*prmax; 947cf1a0672SHong Zhang if (!nlnk_max || nlnk_max > pN) nlnk_max = pN; 948cf1a0672SHong Zhang #if defined(DEBUG_MATMATMULT) 949cf1a0672SHong Zhang if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] pN %d; nlnk_max %d; Armax %d+%d=%d; Prmax Max(%d,%d)=%d\n",rank,pN,nlnk_max,ad->rmax,ao->rmax,armax,p_loc->rmax,p_oth->rmax,prmax); 950cf1a0672SHong Zhang #endif 951cf1a0672SHong Zhang ierr = PetscLLCondensedCreate(nlnk_max,pN,lnk,nlnk,lnkbt); 952cf1a0672SHong Zhang 953cf1a0672SHong Zhang /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */ 954cf1a0672SHong Zhang ierr = PetscFreeSpaceGet((PetscInt)(fill*(adi[am]+aoi[am]+pi_loc[pm])),&free_space);CHKERRQ(ierr); 955cf1a0672SHong Zhang current_space = free_space; 956cf1a0672SHong Zhang 957cf1a0672SHong Zhang ierr = MatPreallocateInitialize(comm,am,pn,dnz,onz);CHKERRQ(ierr); 958cf1a0672SHong Zhang for (i=0; i<am; i++) { 959cf1a0672SHong Zhang apnz = 0; 960cf1a0672SHong Zhang /* diagonal portion of A */ 961cf1a0672SHong Zhang nzi = adi[i+1] - adi[i]; 962cf1a0672SHong Zhang for (j=0; j<nzi; j++){ 963cf1a0672SHong Zhang row = *adj++; 964cf1a0672SHong Zhang pnz = pi_loc[row+1] - pi_loc[row]; 965cf1a0672SHong Zhang Jptr = pj_loc + pi_loc[row]; 966cf1a0672SHong Zhang /* add non-zero cols of P into the sorted linked list lnk */ 967cf1a0672SHong Zhang ierr = PetscLLCondensedAddSorted(nlnk_max,pN,pnz,Jptr,nlnk,lnk,lnkbt);CHKERRQ(ierr); 968cf1a0672SHong Zhang } 969cf1a0672SHong Zhang /* off-diagonal portion of A */ 970cf1a0672SHong Zhang nzi = aoi[i+1] - aoi[i]; 971cf1a0672SHong Zhang for (j=0; j<nzi; j++){ 972cf1a0672SHong Zhang row = *aoj++; 973cf1a0672SHong Zhang pnz = pi_oth[row+1] - pi_oth[row]; 974cf1a0672SHong Zhang Jptr = pj_oth + pi_oth[row]; 975cf1a0672SHong Zhang ierr = PetscLLCondensedAddSorted(nlnk_max,pN,pnz,Jptr,nlnk,lnk,lnkbt);CHKERRQ(ierr); 976cf1a0672SHong Zhang } 977cf1a0672SHong Zhang 978cf1a0672SHong Zhang apnz = *nlnk; 979cf1a0672SHong Zhang api[i+1] = api[i] + apnz; 980cf1a0672SHong Zhang if (apnz > apnz_max) apnz_max = apnz; 981cf1a0672SHong Zhang 982cf1a0672SHong Zhang /* if free space is not available, double the total space in the list */ 983cf1a0672SHong Zhang if (current_space->local_remaining<apnz) { 984cf1a0672SHong Zhang ierr = PetscFreeSpaceGet(apnz+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 985cf1a0672SHong Zhang nspacedouble++; 986cf1a0672SHong Zhang } 987cf1a0672SHong Zhang 988cf1a0672SHong Zhang /* Copy data into free space, then initialize lnk */ 989cf1a0672SHong Zhang ierr = PetscLLCondensedClean(nlnk_max,pN,apnz,current_space->array,nlnk,lnk,lnkbt);CHKERRQ(ierr); 990cf1a0672SHong Zhang ierr = MatPreallocateSet(i+rstart,apnz,current_space->array,dnz,onz);CHKERRQ(ierr); 991cf1a0672SHong Zhang current_space->array += apnz; 992cf1a0672SHong Zhang current_space->local_used += apnz; 993cf1a0672SHong Zhang current_space->local_remaining -= apnz; 994cf1a0672SHong Zhang } 995cf1a0672SHong Zhang 996cf1a0672SHong Zhang /* Allocate space for apj, initialize apj, and */ 997cf1a0672SHong Zhang /* destroy list of free space and other temporary array(s) */ 998cf1a0672SHong Zhang ierr = PetscMalloc((api[am]+1)*sizeof(PetscInt),&ptap->apj);CHKERRQ(ierr); 999cf1a0672SHong Zhang apj = ptap->apj; 1000cf1a0672SHong Zhang ierr = PetscFreeSpaceContiguous(&free_space,ptap->apj);CHKERRQ(ierr); 1001cf1a0672SHong Zhang ierr = PetscLLCondensedDestroy(lnk,lnkbt);CHKERRQ(ierr); 1002cf1a0672SHong Zhang #if defined(DEBUG_MATMATMULT) 1003cf1a0672SHong Zhang if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] AP is done..., apnz_max %d\n",rank,apnz_max); 1004cf1a0672SHong Zhang #endif 1005cf1a0672SHong Zhang 1006cf1a0672SHong Zhang /* create and assemble symbolic parallel matrix Cmpi */ 1007cf1a0672SHong Zhang /*----------------------------------------------------*/ 1008cf1a0672SHong Zhang ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr); 1009cf1a0672SHong Zhang ierr = MatSetSizes(Cmpi,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 1010cf1a0672SHong Zhang ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr); 1011cf1a0672SHong Zhang ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr); 1012cf1a0672SHong Zhang ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 1013cf1a0672SHong Zhang ierr = MatSetBlockSize(Cmpi,1);CHKERRQ(ierr); 1014cf1a0672SHong Zhang 101529825010SHong Zhang /* malloc apa for assembly Cmpi */ 1016cf1a0672SHong Zhang ierr = PetscMalloc(apnz_max*sizeof(PetscScalar),&apa);CHKERRQ(ierr); 1017cf1a0672SHong Zhang ierr = PetscMemzero(apa,apnz_max*sizeof(PetscScalar));CHKERRQ(ierr); 101829825010SHong Zhang ptap->apa = apa; 1019cf1a0672SHong Zhang for (i=0; i<am; i++){ 1020cf1a0672SHong Zhang row = i + rstart; 1021cf1a0672SHong Zhang apnz = api[i+1] - api[i]; 1022cf1a0672SHong Zhang ierr = MatSetValues(Cmpi,1,&row,apnz,apj,apa,INSERT_VALUES);CHKERRQ(ierr); 1023cf1a0672SHong Zhang apj += apnz; 1024cf1a0672SHong Zhang } 1025cf1a0672SHong Zhang ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1026cf1a0672SHong Zhang ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1027cf1a0672SHong Zhang 1028cf1a0672SHong Zhang ptap->destroy = Cmpi->ops->destroy; 1029cf1a0672SHong Zhang ptap->duplicate = Cmpi->ops->duplicate; 1030cf1a0672SHong Zhang Cmpi->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ_Scalable; 1031cf1a0672SHong Zhang Cmpi->ops->destroy = MatDestroy_MPIAIJ_MatMatMult; 1032cf1a0672SHong Zhang Cmpi->ops->duplicate = MatDuplicate_MPIAIJ_MatMatMult; 1033cf1a0672SHong Zhang 1034cf1a0672SHong Zhang /* attach the supporting struct to Cmpi for reuse */ 1035cf1a0672SHong Zhang c = (Mat_MPIAIJ*)Cmpi->data; 1036cf1a0672SHong Zhang c->ptap = ptap; 1037cf1a0672SHong Zhang 1038cf1a0672SHong Zhang *C = Cmpi; 1039cf1a0672SHong Zhang 1040cf1a0672SHong Zhang /* set MatInfo */ 1041cf1a0672SHong Zhang afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]) + 1.e-5; 1042cf1a0672SHong Zhang if (afill < 1.0) afill = 1.0; 1043cf1a0672SHong Zhang Cmpi->info.mallocs = nspacedouble; 1044cf1a0672SHong Zhang Cmpi->info.fill_ratio_given = fill; 1045cf1a0672SHong Zhang Cmpi->info.fill_ratio_needed = afill; 1046cf1a0672SHong Zhang 1047cf1a0672SHong Zhang #if defined(PETSC_USE_INFO) 1048cf1a0672SHong Zhang if (api[am]) { 1049cf1a0672SHong Zhang ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr); 1050cf1a0672SHong Zhang ierr = PetscInfo1(Cmpi,"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr); 1051cf1a0672SHong Zhang } else { 1052cf1a0672SHong Zhang ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr); 1053cf1a0672SHong Zhang } 1054cf1a0672SHong Zhang #endif 10551634b032SHong Zhang PetscFunctionReturn(0); 10561634b032SHong Zhang } 1057