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 */ 14e745005bSHong Zhang /* 151c7d5954SHong Zhang #define DEBUG_MATTrMatMult 16e745005bSHong Zhang */ 171634b032SHong Zhang 182499ec78SHong Zhang #undef __FUNCT__ 192499ec78SHong Zhang #define __FUNCT__ "MatMatMult_MPIAIJ_MPIAIJ" 202499ec78SHong Zhang PetscErrorCode MatMatMult_MPIAIJ_MPIAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill, Mat *C) 212499ec78SHong Zhang { 222499ec78SHong Zhang PetscErrorCode ierr; 232499ec78SHong Zhang 242499ec78SHong Zhang PetscFunctionBegin; 252499ec78SHong Zhang if (scall == MAT_INITIAL_MATRIX){ 26598bc09dSHong Zhang ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr); 27a1a4e84aSHong Zhang ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ(A,B,fill,C);CHKERRQ(ierr); 28598bc09dSHong Zhang ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr); 292499ec78SHong Zhang } 30598bc09dSHong Zhang ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr); 31598bc09dSHong Zhang ierr = (*(*C)->ops->matmultnumeric)(A,B,*C);CHKERRQ(ierr); 32598bc09dSHong Zhang ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr); 332499ec78SHong Zhang PetscFunctionReturn(0); 342499ec78SHong Zhang } 352499ec78SHong Zhang 36f33d1a9aSHong Zhang #undef __FUNCT__ 37776b82aeSLisandro Dalcin #define __FUNCT__ "PetscContainerDestroy_Mat_MatMatMultMPI" 38776b82aeSLisandro Dalcin PetscErrorCode PetscContainerDestroy_Mat_MatMatMultMPI(void *ptr) 39f33d1a9aSHong Zhang { 40f33d1a9aSHong Zhang PetscErrorCode ierr; 419649163dSHong Zhang Mat_MatMatMultMPI *mult=(Mat_MatMatMultMPI*)ptr; 42f33d1a9aSHong Zhang 43f33d1a9aSHong Zhang PetscFunctionBegin; 446bf464f9SBarry Smith ierr = ISDestroy(&mult->isrowa);CHKERRQ(ierr); 456bf464f9SBarry Smith ierr = ISDestroy(&mult->isrowb);CHKERRQ(ierr); 466bf464f9SBarry Smith ierr = ISDestroy(&mult->iscolb);CHKERRQ(ierr); 476bf464f9SBarry Smith ierr = MatDestroy(&mult->C_seq);CHKERRQ(ierr); 486bf464f9SBarry Smith ierr = MatDestroy(&mult->A_loc);CHKERRQ(ierr); 496bf464f9SBarry Smith ierr = MatDestroy(&mult->B_seq);CHKERRQ(ierr); 50f33d1a9aSHong Zhang ierr = PetscFree(mult);CHKERRQ(ierr); 51f33d1a9aSHong Zhang PetscFunctionReturn(0); 52f33d1a9aSHong Zhang } 53f33d1a9aSHong Zhang 54598bc09dSHong Zhang #undef __FUNCT__ 55a1a4e84aSHong Zhang #define __FUNCT__ "MatDestroy_MPIAIJ_MatMatMult" 56a1a4e84aSHong Zhang PetscErrorCode MatDestroy_MPIAIJ_MatMatMult(Mat A) 57598bc09dSHong Zhang { 58598bc09dSHong Zhang PetscErrorCode ierr; 59598bc09dSHong Zhang Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data; 60598bc09dSHong Zhang Mat_PtAPMPI *ptap=a->ptap; 61598bc09dSHong Zhang 62598bc09dSHong Zhang PetscFunctionBegin; 63b7f45c76SHong Zhang ierr = PetscFree2(ptap->startsj_s,ptap->startsj_r);CHKERRQ(ierr); 64598bc09dSHong Zhang ierr = PetscFree(ptap->bufa);CHKERRQ(ierr); 65a1a4e84aSHong Zhang ierr = MatDestroy(&ptap->P_loc);CHKERRQ(ierr); 66a1a4e84aSHong Zhang ierr = MatDestroy(&ptap->P_oth);CHKERRQ(ierr); 67a1a4e84aSHong Zhang ierr = PetscFree(ptap->api);CHKERRQ(ierr); 68a1a4e84aSHong Zhang ierr = PetscFree(ptap->apj);CHKERRQ(ierr); 69598bc09dSHong Zhang ierr = PetscFree(ptap->apa);CHKERRQ(ierr); 70598bc09dSHong Zhang ierr = ptap->destroy(A);CHKERRQ(ierr); 71598bc09dSHong Zhang ierr = PetscFree(ptap);CHKERRQ(ierr); 72598bc09dSHong Zhang PetscFunctionReturn(0); 73598bc09dSHong Zhang } 74598bc09dSHong Zhang 752499ec78SHong Zhang #undef __FUNCT__ 76a1a4e84aSHong Zhang #define __FUNCT__ "MatDestroy_MPIAIJ_MatMatMult_32" 77a1a4e84aSHong Zhang PetscErrorCode MatDestroy_MPIAIJ_MatMatMult_32(Mat A) 782499ec78SHong Zhang { 792499ec78SHong Zhang PetscErrorCode ierr; 80776b82aeSLisandro Dalcin PetscContainer container; 819649163dSHong Zhang Mat_MatMatMultMPI *mult=PETSC_NULL; 822499ec78SHong Zhang 832499ec78SHong Zhang PetscFunctionBegin; 849649163dSHong Zhang ierr = PetscObjectQuery((PetscObject)A,"Mat_MatMatMultMPI",(PetscObject *)&container);CHKERRQ(ierr); 8529825010SHong Zhang if (!container) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_PLIB,"Container does not exist"); 86776b82aeSLisandro Dalcin ierr = PetscContainerGetPointer(container,(void **)&mult);CHKERRQ(ierr); 87dce485f0SBarry Smith A->ops->destroy = mult->destroy; 88bf0cc555SLisandro Dalcin A->ops->duplicate = mult->duplicate; 89bf0cc555SLisandro Dalcin if (A->ops->destroy) { 90992edd73SBarry Smith ierr = (*A->ops->destroy)(A);CHKERRQ(ierr); 91bf0cc555SLisandro Dalcin } 92bf0cc555SLisandro Dalcin ierr = PetscObjectCompose((PetscObject)A,"Mat_MatMatMultMPI",0);CHKERRQ(ierr); 932499ec78SHong Zhang PetscFunctionReturn(0); 942499ec78SHong Zhang } 952499ec78SHong Zhang 962499ec78SHong Zhang #undef __FUNCT__ 97a1a4e84aSHong Zhang #define __FUNCT__ "MatDuplicate_MPIAIJ_MatMatMult_32" 98a1a4e84aSHong Zhang PetscErrorCode MatDuplicate_MPIAIJ_MatMatMult_32(Mat A, MatDuplicateOption op, Mat *M) 99b9c92825SBarry Smith { 100abf897f8SHong Zhang PetscErrorCode ierr; 101abf897f8SHong Zhang Mat_MatMatMultMPI *mult; 102776b82aeSLisandro Dalcin PetscContainer container; 103abf897f8SHong Zhang 104abf897f8SHong Zhang PetscFunctionBegin; 105abf897f8SHong Zhang ierr = PetscObjectQuery((PetscObject)A,"Mat_MatMatMultMPI",(PetscObject *)&container);CHKERRQ(ierr); 10629825010SHong Zhang if (!container) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_PLIB,"Container does not exist"); 107776b82aeSLisandro Dalcin ierr = PetscContainerGetPointer(container,(void **)&mult);CHKERRQ(ierr); 1085c088420SHong Zhang /* Note: the container is not duplicated, because it requires deep copying of 1095c088420SHong Zhang several large data sets (see PetscContainerDestroy_Mat_MatMatMultMPI()). 1105c088420SHong Zhang These data sets are only used for repeated calling of MatMatMultNumeric(). 1115c088420SHong Zhang *M is unlikely being used in this way. Thus we create *M with pure mpiaij format */ 112dce485f0SBarry Smith ierr = (*mult->duplicate)(A,op,M);CHKERRQ(ierr); 113dce485f0SBarry Smith (*M)->ops->destroy = mult->destroy; /* = MatDestroy_MPIAIJ, *M doesn't duplicate A's container! */ 114dce485f0SBarry Smith (*M)->ops->duplicate = mult->duplicate; /* = MatDuplicate_MPIAIJ */ 115abf897f8SHong Zhang PetscFunctionReturn(0); 116abf897f8SHong Zhang } 117abf897f8SHong Zhang 118abf897f8SHong Zhang #undef __FUNCT__ 119a1a4e84aSHong Zhang #define __FUNCT__ "MatDuplicate_MPIAIJ_MatMatMult" 120a1a4e84aSHong Zhang PetscErrorCode MatDuplicate_MPIAIJ_MatMatMult(Mat A, MatDuplicateOption op, Mat *M) 1214ae0847bSHong Zhang { 1224ae0847bSHong Zhang PetscErrorCode ierr; 1234ae0847bSHong Zhang Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data; 1244ae0847bSHong Zhang Mat_PtAPMPI *ptap=a->ptap; 1254ae0847bSHong Zhang 1264ae0847bSHong Zhang PetscFunctionBegin; 1274ae0847bSHong Zhang ierr = (*ptap->duplicate)(A,op,M);CHKERRQ(ierr); 1284ae0847bSHong Zhang (*M)->ops->destroy = ptap->destroy; /* = MatDestroy_MPIAIJ, *M doesn't duplicate A's special structure! */ 1294ae0847bSHong Zhang (*M)->ops->duplicate = ptap->duplicate; /* = MatDuplicate_MPIAIJ */ 1304ae0847bSHong Zhang PetscFunctionReturn(0); 1314ae0847bSHong Zhang } 1324ae0847bSHong Zhang 1334ae0847bSHong Zhang #undef __FUNCT__ 134a1a4e84aSHong Zhang #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIAIJ" 135a1a4e84aSHong Zhang PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C) 136598bc09dSHong Zhang { 137598bc09dSHong Zhang PetscErrorCode ierr; 1384ae0847bSHong Zhang Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data,*c=(Mat_MPIAIJ*)C->data; 139598bc09dSHong Zhang Mat_SeqAIJ *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data; 140598bc09dSHong Zhang Mat_SeqAIJ *cd=(Mat_SeqAIJ*)(c->A)->data,*co=(Mat_SeqAIJ*)(c->B)->data; 141598bc09dSHong Zhang PetscInt *adi=ad->i,*adj,*aoi=ao->i,*aoj; 142598bc09dSHong Zhang PetscScalar *ada,*aoa,*cda=cd->a,*coa=co->a; 143598bc09dSHong Zhang Mat_SeqAIJ *p_loc,*p_oth; 144598bc09dSHong Zhang PetscInt *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pj; 145598bc09dSHong Zhang PetscScalar *pa_loc,*pa_oth,*pa,*apa,valtmp,*ca; 146598bc09dSHong Zhang PetscInt cm=C->rmap->n,anz,pnz; 147598bc09dSHong Zhang Mat_PtAPMPI *ptap=c->ptap; 14865a57a32SSean Farley PetscInt *api,*apj,*apJ,i,j,k,row; 14929825010SHong Zhang PetscInt cstart=C->cmap->rstart; 150598bc09dSHong Zhang PetscInt cdnz,conz,k0,k1; 1517c2e2f26SHong Zhang #if defined(DEBUG_MATMATMULT) 1527c2e2f26SHong Zhang PetscMPIInt rank=a->rank; 1537c2e2f26SHong Zhang #endif 154598bc09dSHong Zhang 155598bc09dSHong Zhang PetscFunctionBegin; 1567c2e2f26SHong Zhang #if defined(DEBUG_MATMATMULT) 1577c2e2f26SHong Zhang if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] call MatMatMultNumeric_MPIAIJ_MPIAIJ()...\n",rank); 1587c2e2f26SHong Zhang #endif 1597c2e2f26SHong Zhang 160a1a4e84aSHong Zhang /* 1) get P_oth = ptap->P_oth and P_loc = ptap->P_loc */ 161598bc09dSHong Zhang /*-----------------------------------------------------*/ 162a1a4e84aSHong Zhang /* update numerical values of P_oth and P_loc */ 163b7f45c76SHong Zhang ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 164a1a4e84aSHong Zhang ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 1656291f5a6SHong Zhang #if defined(DEBUG_MATMATMULT) 1666291f5a6SHong Zhang if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] got P_oth and P_loc...\n",rank); 1676291f5a6SHong Zhang #endif 168598bc09dSHong Zhang 169598bc09dSHong Zhang /* 2) compute numeric C_loc = A_loc*P = Ad*P_loc + Ao*P_oth */ 170598bc09dSHong Zhang /*----------------------------------------------------------*/ 171598bc09dSHong Zhang /* get data from symbolic products */ 172a1a4e84aSHong Zhang p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 173a1a4e84aSHong Zhang p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 174598bc09dSHong Zhang pi_loc=p_loc->i; pj_loc=p_loc->j; pa_loc=p_loc->a; 175598bc09dSHong Zhang pi_oth=p_oth->i; pj_oth=p_oth->j; pa_oth=p_oth->a; 176598bc09dSHong Zhang 177598bc09dSHong Zhang /* get apa for storing dense row A[i,:]*P */ 178598bc09dSHong Zhang apa = ptap->apa; 179598bc09dSHong Zhang 18029825010SHong Zhang api = ptap->api; 18129825010SHong Zhang apj = ptap->apj; 182598bc09dSHong Zhang for (i=0; i<cm; i++) { 183598bc09dSHong Zhang /* diagonal portion of A */ 184598bc09dSHong Zhang anz = adi[i+1] - adi[i]; 185598bc09dSHong Zhang adj = ad->j + adi[i]; 186598bc09dSHong Zhang ada = ad->a + adi[i]; 187598bc09dSHong Zhang for (j=0; j<anz; j++) { 188598bc09dSHong Zhang row = adj[j]; 189598bc09dSHong Zhang pnz = pi_loc[row+1] - pi_loc[row]; 190598bc09dSHong Zhang pj = pj_loc + pi_loc[row]; 191598bc09dSHong Zhang pa = pa_loc + pi_loc[row]; 192598bc09dSHong Zhang 193598bc09dSHong Zhang /* perform dense axpy */ 194598bc09dSHong Zhang valtmp = ada[j]; 195598bc09dSHong Zhang for (k=0; k<pnz; k++){ 196598bc09dSHong Zhang apa[pj[k]] += valtmp*pa[k]; 197598bc09dSHong Zhang } 198598bc09dSHong Zhang ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr); 199598bc09dSHong Zhang } 200598bc09dSHong Zhang 201598bc09dSHong Zhang /* off-diagonal portion of A */ 202598bc09dSHong Zhang anz = aoi[i+1] - aoi[i]; 203598bc09dSHong Zhang aoj = ao->j + aoi[i]; 204598bc09dSHong Zhang aoa = ao->a + aoi[i]; 205598bc09dSHong Zhang for (j=0; j<anz; j++) { 206598bc09dSHong Zhang row = aoj[j]; 207598bc09dSHong Zhang pnz = pi_oth[row+1] - pi_oth[row]; 208598bc09dSHong Zhang pj = pj_oth + pi_oth[row]; 209598bc09dSHong Zhang pa = pa_oth + pi_oth[row]; 210598bc09dSHong Zhang 211598bc09dSHong Zhang /* perform dense axpy */ 212598bc09dSHong Zhang valtmp = aoa[j]; 213598bc09dSHong Zhang for (k=0; k<pnz; k++){ 214598bc09dSHong Zhang apa[pj[k]] += valtmp*pa[k]; 215598bc09dSHong Zhang } 216598bc09dSHong Zhang ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr); 217598bc09dSHong Zhang } 218598bc09dSHong Zhang 219598bc09dSHong Zhang /* set values in C */ 220598bc09dSHong Zhang apJ = apj + api[i]; 221598bc09dSHong Zhang cdnz = cd->i[i+1] - cd->i[i]; 222598bc09dSHong Zhang conz = co->i[i+1] - co->i[i]; 223598bc09dSHong Zhang 224598bc09dSHong Zhang /* 1st off-diagoanl part of C */ 225598bc09dSHong Zhang ca = coa + co->i[i]; 226598bc09dSHong Zhang k = 0; 227598bc09dSHong Zhang for (k0=0; k0<conz; k0++){ 228598bc09dSHong Zhang if (apJ[k] >= cstart) break; 229598bc09dSHong Zhang ca[k0] = apa[apJ[k]]; 230598bc09dSHong Zhang apa[apJ[k]] = 0.0; 231598bc09dSHong Zhang k++; 232598bc09dSHong Zhang } 233598bc09dSHong Zhang 234598bc09dSHong Zhang /* diagonal part of C */ 235598bc09dSHong Zhang ca = cda + cd->i[i]; 236598bc09dSHong Zhang for (k1=0; k1<cdnz; k1++){ 237598bc09dSHong Zhang ca[k1] = apa[apJ[k]]; 238598bc09dSHong Zhang apa[apJ[k]] = 0.0; 239598bc09dSHong Zhang k++; 240598bc09dSHong Zhang } 241598bc09dSHong Zhang 242598bc09dSHong Zhang /* 2nd off-diagoanl part of C */ 243598bc09dSHong Zhang ca = coa + co->i[i]; 244598bc09dSHong Zhang for (; k0<conz; k0++){ 245598bc09dSHong Zhang ca[k0] = apa[apJ[k]]; 246598bc09dSHong Zhang apa[apJ[k]] = 0.0; 247598bc09dSHong Zhang k++; 248598bc09dSHong Zhang } 249598bc09dSHong Zhang } 250598bc09dSHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 251598bc09dSHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 252598bc09dSHong Zhang PetscFunctionReturn(0); 253598bc09dSHong Zhang } 254598bc09dSHong Zhang 255598bc09dSHong Zhang #undef __FUNCT__ 256a1a4e84aSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ" 257a1a4e84aSHong Zhang PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C) 258598bc09dSHong Zhang { 259598bc09dSHong Zhang PetscErrorCode ierr; 260598bc09dSHong Zhang MPI_Comm comm=((PetscObject)A)->comm; 261bfcd1a73SHong Zhang Mat Cmpi; 262598bc09dSHong Zhang Mat_PtAPMPI *ptap; 263598bc09dSHong Zhang PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 2644ae0847bSHong Zhang Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data,*c; 265bfcd1a73SHong Zhang Mat_SeqAIJ *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_loc,*p_oth; 2664ae0847bSHong Zhang PetscInt *pi_loc,*pj_loc,*pi_oth,*pj_oth,*dnz,*onz; 2674ae0847bSHong Zhang PetscInt *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,rstart=A->rmap->rstart; 268d14fa243SHong Zhang PetscInt *lnk,i,pnz,row,*api,*apj,*Jptr,apnz,nspacedouble=0,j,nzi; 269bfcd1a73SHong Zhang PetscInt am=A->rmap->n,pN=P->cmap->N,pn=P->cmap->n,pm=P->rmap->n; 270598bc09dSHong Zhang PetscBT lnkbt; 271598bc09dSHong Zhang PetscScalar *apa; 272bfcd1a73SHong Zhang PetscReal afill; 273cf1a0672SHong Zhang PetscBool scalable=PETSC_FALSE; 274f84c536bSHong Zhang PetscInt nlnk_max,armax,prmax; 2757c2e2f26SHong Zhang #if defined(DEBUG_MATMATMULT) 2767c2e2f26SHong Zhang PetscMPIInt rank=a->rank; 2777c2e2f26SHong Zhang #endif 278598bc09dSHong Zhang 279598bc09dSHong Zhang PetscFunctionBegin; 280a1a4e84aSHong Zhang if (A->cmap->rstart != P->rmap->rstart || A->cmap->rend != P->rmap->rend){ 281cf1a0672SHong 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); 282a1a4e84aSHong Zhang } 283152983bfSHong Zhang 284152983bfSHong Zhang ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr); 285cf1a0672SHong Zhang 286cf1a0672SHong Zhang ierr = PetscOptionsBool("-matmatmult_32","Use a scalable but slower C=A*B","",scalable,&scalable,PETSC_NULL);CHKERRQ(ierr); 287cf1a0672SHong Zhang if (scalable){ 288cf1a0672SHong Zhang ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable_32(A,P,fill,C);;CHKERRQ(ierr); 289a1a4e84aSHong Zhang PetscFunctionReturn(0); 290a1a4e84aSHong Zhang } 291cf1a0672SHong Zhang ierr = PetscOptionsBool("-matmatmult_scalable","Use a scalable but slower C=A*B","",scalable,&scalable,PETSC_NULL);CHKERRQ(ierr); 292cf1a0672SHong Zhang if (scalable){ 29325023028SHong Zhang ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable(A,P,fill,C);;CHKERRQ(ierr); 2941634b032SHong Zhang PetscFunctionReturn(0); 2951634b032SHong Zhang } 296152983bfSHong Zhang ierr = PetscOptionsEnd();CHKERRQ(ierr); 297a1a4e84aSHong Zhang 2987c2e2f26SHong Zhang #if defined(DEBUG_MATMATMULT) 2997c2e2f26SHong Zhang if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] call MatMatMultSymbolic_MPIAIJ_MPIAIJ()...\n",rank); 3007c2e2f26SHong Zhang #endif 3017c2e2f26SHong Zhang 302598bc09dSHong Zhang /* create struct Mat_PtAPMPI and attached it to C later */ 303598bc09dSHong Zhang ierr = PetscNew(Mat_PtAPMPI,&ptap);CHKERRQ(ierr); 304598bc09dSHong Zhang 305598bc09dSHong Zhang /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */ 306b7f45c76SHong Zhang ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 3076291f5a6SHong Zhang #if defined(DEBUG_MATMATMULT) 3086291f5a6SHong Zhang if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] P_oth is done...\n",rank); 3096291f5a6SHong Zhang #endif 310598bc09dSHong Zhang /* get P_loc by taking all local rows of P */ 311a1a4e84aSHong Zhang ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 312598bc09dSHong Zhang 313a1a4e84aSHong Zhang p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 314a1a4e84aSHong Zhang p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 315598bc09dSHong Zhang pi_loc = p_loc->i; pj_loc = p_loc->j; 316598bc09dSHong Zhang pi_oth = p_oth->i; pj_oth = p_oth->j; 317598bc09dSHong Zhang 3186291f5a6SHong Zhang #if defined(DEBUG_MATMATMULT) 319e9e4536cSHong 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); 3206291f5a6SHong Zhang #endif 3216291f5a6SHong Zhang 322598bc09dSHong Zhang /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */ 323598bc09dSHong Zhang /*-------------------------------------------------------------------*/ 324598bc09dSHong Zhang ierr = PetscMalloc((am+2)*sizeof(PetscInt),&api);CHKERRQ(ierr); 325a1a4e84aSHong Zhang ptap->api = api; 326598bc09dSHong Zhang api[0] = 0; 327598bc09dSHong Zhang 328598bc09dSHong Zhang /* create and initialize a linked list */ 329f84c536bSHong Zhang armax = ad->rmax+ao->rmax; 330f84c536bSHong Zhang prmax = PetscMax(p_loc->rmax,p_oth->rmax); 331f84c536bSHong Zhang nlnk_max = armax*prmax; 332f84c536bSHong Zhang if (!nlnk_max || nlnk_max > pN) nlnk_max = pN; 333f84c536bSHong Zhang #if defined(DEBUG_MATMATMULT) 334f84c536bSHong 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); 335f84c536bSHong Zhang #endif 3360d3134e9SHong Zhang ierr = PetscLLCondensedCreate(nlnk_max,pN,&lnk,&lnkbt);CHKERRQ(ierr); 337598bc09dSHong Zhang 338bfcd1a73SHong Zhang /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */ 339bfcd1a73SHong Zhang ierr = PetscFreeSpaceGet((PetscInt)(fill*(adi[am]+aoi[am]+pi_loc[pm])),&free_space);CHKERRQ(ierr); 340598bc09dSHong Zhang current_space = free_space; 341598bc09dSHong Zhang 342598bc09dSHong Zhang ierr = MatPreallocateInitialize(comm,am,pn,dnz,onz);CHKERRQ(ierr); 343598bc09dSHong Zhang for (i=0; i<am; i++) { 344598bc09dSHong Zhang apnz = 0; 345598bc09dSHong Zhang /* diagonal portion of A */ 346598bc09dSHong Zhang nzi = adi[i+1] - adi[i]; 347598bc09dSHong Zhang for (j=0; j<nzi; j++){ 348598bc09dSHong Zhang row = *adj++; 349598bc09dSHong Zhang pnz = pi_loc[row+1] - pi_loc[row]; 350598bc09dSHong Zhang Jptr = pj_loc + pi_loc[row]; 351598bc09dSHong Zhang /* add non-zero cols of P into the sorted linked list lnk */ 3521fbbb359SHong Zhang ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 353598bc09dSHong Zhang } 354598bc09dSHong Zhang /* off-diagonal portion of A */ 355598bc09dSHong Zhang nzi = aoi[i+1] - aoi[i]; 356598bc09dSHong Zhang for (j=0; j<nzi; j++){ 357598bc09dSHong Zhang row = *aoj++; 358598bc09dSHong Zhang pnz = pi_oth[row+1] - pi_oth[row]; 359598bc09dSHong Zhang Jptr = pj_oth + pi_oth[row]; 3601fbbb359SHong Zhang ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 361598bc09dSHong Zhang } 362598bc09dSHong Zhang 363d14fa243SHong Zhang apnz = lnk[0]; 364598bc09dSHong Zhang api[i+1] = api[i] + apnz; 365598bc09dSHong Zhang 366598bc09dSHong Zhang /* if free space is not available, double the total space in the list */ 367598bc09dSHong Zhang if (current_space->local_remaining<apnz) { 368598bc09dSHong Zhang ierr = PetscFreeSpaceGet(apnz+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 369598bc09dSHong Zhang nspacedouble++; 370598bc09dSHong Zhang } 371598bc09dSHong Zhang 372598bc09dSHong Zhang /* Copy data into free space, then initialize lnk */ 3731fbbb359SHong Zhang ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 374598bc09dSHong Zhang ierr = MatPreallocateSet(i+rstart,apnz,current_space->array,dnz,onz);CHKERRQ(ierr); 375598bc09dSHong Zhang current_space->array += apnz; 376598bc09dSHong Zhang current_space->local_used += apnz; 377598bc09dSHong Zhang current_space->local_remaining -= apnz; 378598bc09dSHong Zhang } 379598bc09dSHong Zhang 380598bc09dSHong Zhang /* Allocate space for apj, initialize apj, and */ 381598bc09dSHong Zhang /* destroy list of free space and other temporary array(s) */ 382a1a4e84aSHong Zhang ierr = PetscMalloc((api[am]+1)*sizeof(PetscInt),&ptap->apj);CHKERRQ(ierr); 383a1a4e84aSHong Zhang apj = ptap->apj; 384a1a4e84aSHong Zhang ierr = PetscFreeSpaceContiguous(&free_space,ptap->apj);CHKERRQ(ierr); 385598bc09dSHong Zhang ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 3866291f5a6SHong Zhang #if defined(DEBUG_MATMATMULT) 3876291f5a6SHong Zhang if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] AP is done...\n",rank); 3886291f5a6SHong Zhang #endif 389598bc09dSHong Zhang 390f84c536bSHong Zhang /* malloc apa to store dense row A[i,:]*P */ 391f84c536bSHong Zhang ierr = PetscMalloc(pN*sizeof(PetscScalar),&apa);CHKERRQ(ierr); 392f84c536bSHong Zhang ierr = PetscMemzero(apa,pN*sizeof(PetscScalar));CHKERRQ(ierr); 393f84c536bSHong Zhang ptap->apa = apa; 394f84c536bSHong Zhang #if defined(DEBUG_MATMATMULT) 395f84c536bSHong Zhang if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] Malloc apa pN %D is done...\n",rank,pN); 396f84c536bSHong Zhang #endif 397f84c536bSHong Zhang 398bfcd1a73SHong Zhang /* create and assemble symbolic parallel matrix Cmpi */ 399598bc09dSHong Zhang /*----------------------------------------------------*/ 400bfcd1a73SHong Zhang ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr); 401bfcd1a73SHong Zhang ierr = MatSetSizes(Cmpi,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 402bfcd1a73SHong Zhang ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr); 403bfcd1a73SHong Zhang ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr); 404598bc09dSHong Zhang ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 405598bc09dSHong Zhang for (i=0; i<am; i++){ 406598bc09dSHong Zhang row = i + rstart; 407598bc09dSHong Zhang apnz = api[i+1] - api[i]; 408bfcd1a73SHong Zhang ierr = MatSetValues(Cmpi,1,&row,apnz,apj,apa,INSERT_VALUES);CHKERRQ(ierr); 409598bc09dSHong Zhang apj += apnz; 410598bc09dSHong Zhang } 411bfcd1a73SHong Zhang ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 412bfcd1a73SHong Zhang ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 413598bc09dSHong Zhang 414bfcd1a73SHong Zhang ptap->destroy = Cmpi->ops->destroy; 415bfcd1a73SHong Zhang ptap->duplicate = Cmpi->ops->duplicate; 416bfcd1a73SHong Zhang Cmpi->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ; 417bfcd1a73SHong Zhang Cmpi->ops->destroy = MatDestroy_MPIAIJ_MatMatMult; 418bfcd1a73SHong Zhang Cmpi->ops->duplicate = MatDuplicate_MPIAIJ_MatMatMult; 419598bc09dSHong Zhang 420bfcd1a73SHong Zhang /* attach the supporting struct to Cmpi for reuse */ 421bfcd1a73SHong Zhang c = (Mat_MPIAIJ*)Cmpi->data; 422598bc09dSHong Zhang c->ptap = ptap; 423598bc09dSHong Zhang 424bfcd1a73SHong Zhang *C = Cmpi; 425bfcd1a73SHong Zhang 426bfcd1a73SHong Zhang /* set MatInfo */ 427bfcd1a73SHong Zhang afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]) + 1.e-5; 428bfcd1a73SHong Zhang if (afill < 1.0) afill = 1.0; 429bfcd1a73SHong Zhang Cmpi->info.mallocs = nspacedouble; 430bfcd1a73SHong Zhang Cmpi->info.fill_ratio_given = fill; 431bfcd1a73SHong Zhang Cmpi->info.fill_ratio_needed = afill; 432bfcd1a73SHong Zhang 433bfcd1a73SHong Zhang #if defined(PETSC_USE_INFO) 434bfcd1a73SHong Zhang if (api[am]) { 435bfcd1a73SHong Zhang ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr); 436bfcd1a73SHong Zhang ierr = PetscInfo1(Cmpi,"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr); 437bfcd1a73SHong Zhang } else { 438bfcd1a73SHong Zhang ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr); 439bfcd1a73SHong Zhang } 440bfcd1a73SHong Zhang #endif 441598bc09dSHong Zhang PetscFunctionReturn(0); 442598bc09dSHong Zhang } 443598bc09dSHong Zhang 444a1a4e84aSHong Zhang /* implementation used in PETSc-3.2 */ 445a1a4e84aSHong Zhang /* This routine is called ONLY in the case of reusing previously computed symbolic C */ 446598bc09dSHong Zhang #undef __FUNCT__ 447cf1a0672SHong Zhang #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIAIJ_Scalable_32" 448cf1a0672SHong Zhang PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ_Scalable_32(Mat A,Mat B,Mat C) 449a1a4e84aSHong Zhang { 450a1a4e84aSHong Zhang PetscErrorCode ierr; 451a1a4e84aSHong Zhang Mat *seq; 452a1a4e84aSHong Zhang Mat_MatMatMultMPI *mult; 453a1a4e84aSHong Zhang PetscContainer container; 454a1a4e84aSHong Zhang 455a1a4e84aSHong Zhang PetscFunctionBegin; 45615f4e0f4SHong Zhang ierr = PetscObjectQuery((PetscObject)C,"Mat_MatMatMultMPI",(PetscObject *)&container);CHKERRQ(ierr); 45729825010SHong Zhang if (!container) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_PLIB,"Container does not exist"); 45815f4e0f4SHong Zhang ierr = PetscContainerGetPointer(container,(void **)&mult);CHKERRQ(ierr); 45915f4e0f4SHong Zhang 46015f4e0f4SHong Zhang if (mult->skipNumeric){ /* first numeric product is done during symbolic product */ 46115f4e0f4SHong Zhang mult->skipNumeric = PETSC_FALSE; 46215f4e0f4SHong Zhang PetscFunctionReturn(0); 46315f4e0f4SHong Zhang } 4647c2e2f26SHong Zhang #if defined(DEBUG_MATMATMULT) 4657c2e2f26SHong Zhang PetscMPIInt rank; 466fe615159SHong Zhang MPI_Comm comm = ((PetscObject)C)->comm; 467fe615159SHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 468fe615159SHong Zhang ierr = MPI_Barrier(comm);CHKERRQ(ierr); 469cf1a0672SHong Zhang if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] call MatMatMultNumeric_MPIAIJ_MPIAIJ_Scalable_32()...\n",rank); 4707c2e2f26SHong Zhang #endif 47115f4e0f4SHong Zhang 472a1a4e84aSHong Zhang seq = &mult->B_seq; 473a1a4e84aSHong Zhang ierr = MatGetSubMatrices(B,1,&mult->isrowb,&mult->iscolb,MAT_REUSE_MATRIX,&seq);CHKERRQ(ierr); 474a1a4e84aSHong Zhang mult->B_seq = *seq; 475a1a4e84aSHong Zhang 476a1a4e84aSHong Zhang seq = &mult->A_loc; 477a1a4e84aSHong Zhang ierr = MatGetSubMatrices(A,1,&mult->isrowa,&mult->isrowb,MAT_REUSE_MATRIX,&seq);CHKERRQ(ierr); 478a1a4e84aSHong Zhang mult->A_loc = *seq; 479a1a4e84aSHong Zhang 48025023028SHong Zhang ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(mult->A_loc,mult->B_seq,mult->C_seq);CHKERRQ(ierr); 481a1a4e84aSHong Zhang 482a1a4e84aSHong Zhang ierr = PetscObjectReference((PetscObject)mult->C_seq);CHKERRQ(ierr); 483a1a4e84aSHong Zhang ierr = MatMerge(((PetscObject)A)->comm,mult->C_seq,B->cmap->n,MAT_REUSE_MATRIX,&C);CHKERRQ(ierr); 484a1a4e84aSHong Zhang PetscFunctionReturn(0); 485a1a4e84aSHong Zhang } 486a1a4e84aSHong Zhang 4877c2e2f26SHong Zhang /* numeric product is computed as well */ 488a1a4e84aSHong Zhang #undef __FUNCT__ 489cf1a0672SHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable_32" 490cf1a0672SHong Zhang PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable_32(Mat A,Mat B,PetscReal fill,Mat *C) 4912499ec78SHong Zhang { 4922499ec78SHong Zhang PetscErrorCode ierr; 4932499ec78SHong Zhang Mat_MatMatMultMPI *mult; 494776b82aeSLisandro Dalcin PetscContainer container; 495d5821860SHong Zhang Mat AB,*seq; 496d5821860SHong Zhang Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data; 497d5821860SHong Zhang PetscInt *idx,i,start,ncols,nzA,nzB,*cmap,imark; 4987c2e2f26SHong Zhang #if defined(DEBUG_MATMATMULT) 499e9e4536cSHong Zhang MPI_Comm comm = ((PetscObject)A)->comm; 5007c2e2f26SHong Zhang PetscMPIInt rank=a->rank; 5017c2e2f26SHong Zhang #endif 5022499ec78SHong Zhang 5032499ec78SHong Zhang PetscFunctionBegin; 5047c2e2f26SHong Zhang #if defined(DEBUG_MATMATMULT) 505cf1a0672SHong Zhang if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] call MatMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable_32()...\n",rank); 5067c2e2f26SHong Zhang #endif 507d0f46423SBarry Smith if (A->cmap->rstart != B->rmap->rstart || A->cmap->rend != B->rmap->rend){ 508cf1a0672SHong 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); 5092499ec78SHong Zhang } 510598bc09dSHong Zhang 5112499ec78SHong Zhang ierr = PetscNew(Mat_MatMatMultMPI,&mult);CHKERRQ(ierr); 5122499ec78SHong Zhang 513d5821860SHong Zhang /* get isrowb: nonzero col of A */ 514d5821860SHong Zhang start = A->cmap->rstart; 515d5821860SHong Zhang cmap = a->garray; 516d5821860SHong Zhang nzA = a->A->cmap->n; 517d5821860SHong Zhang nzB = a->B->cmap->n; 518d5821860SHong Zhang ierr = PetscMalloc((nzA+nzB)*sizeof(PetscInt), &idx);CHKERRQ(ierr); 519d5821860SHong Zhang ncols = 0; 520d5821860SHong Zhang for (i=0; i<nzB; i++) { /* row < local row index */ 521d5821860SHong Zhang if (cmap[i] < start) idx[ncols++] = cmap[i]; 522d5821860SHong Zhang else break; 523d5821860SHong Zhang } 524d5821860SHong Zhang imark = i; 525d5821860SHong Zhang for (i=0; i<nzA; i++) idx[ncols++] = start + i; /* local rows */ 526d5821860SHong Zhang for (i=imark; i<nzB; i++) idx[ncols++] = cmap[i]; /* row > local row index */ 527d5821860SHong Zhang ierr = ISCreateGeneral(PETSC_COMM_SELF,ncols,idx,PETSC_OWN_POINTER,&mult->isrowb);CHKERRQ(ierr); 528d5821860SHong Zhang ierr = ISCreateStride(PETSC_COMM_SELF,B->cmap->N,0,1,&mult->iscolb);CHKERRQ(ierr); 529d5821860SHong Zhang 530d5821860SHong Zhang /* get isrowa: all local rows of A */ 531d5821860SHong Zhang ierr = ISCreateStride(PETSC_COMM_SELF,A->rmap->n,A->rmap->rstart,1,&mult->isrowa);CHKERRQ(ierr); 532d5821860SHong Zhang 533d5821860SHong Zhang /* Below should go to MatMatMultNumeric_MPIAIJ_MPIAIJ() - How to generate C there? */ 5342499ec78SHong Zhang /* create a seq matrix B_seq = submatrix of B by taking rows of B that equal to nonzero col of A */ 535d5821860SHong Zhang ierr = MatGetSubMatrices(B,1,&mult->isrowb,&mult->iscolb,MAT_INITIAL_MATRIX,&seq);CHKERRQ(ierr); 536d5821860SHong Zhang mult->B_seq = *seq; 537d5821860SHong Zhang ierr = PetscFree(seq);CHKERRQ(ierr); 5386291f5a6SHong Zhang #if defined(DEBUG_MATMATMULT) 5396291f5a6SHong Zhang if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] B_seq is done...\n",rank); 5406291f5a6SHong Zhang #endif 5412499ec78SHong Zhang 5422499ec78SHong Zhang /* create a seq matrix A_seq = submatrix of A by taking all local rows of A */ 543d5821860SHong Zhang ierr = MatGetSubMatrices(A,1,&mult->isrowa,&mult->isrowb,MAT_INITIAL_MATRIX,&seq);CHKERRQ(ierr); 544d5821860SHong Zhang mult->A_loc = *seq; 545d5821860SHong Zhang ierr = PetscFree(seq);CHKERRQ(ierr); 5466291f5a6SHong Zhang #if defined(DEBUG_MATMATMULT) 547e9e4536cSHong Zhang if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] A_loc is done...\n",rank); 5486291f5a6SHong Zhang #endif 5492499ec78SHong Zhang 5502499ec78SHong Zhang /* compute C_seq = A_seq * B_seq */ 55125023028SHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(mult->A_loc,mult->B_seq,fill,&mult->C_seq);CHKERRQ(ierr); 5526291f5a6SHong Zhang #if defined(DEBUG_MATMATMULT) 553e9e4536cSHong Zhang ierr = MPI_Barrier(comm);CHKERRQ(ierr); 554e9e4536cSHong Zhang if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] C_seq Symbolic is done...\n",rank); 555e9e4536cSHong Zhang #endif 55625023028SHong Zhang ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(mult->A_loc,mult->B_seq,mult->C_seq);CHKERRQ(ierr); 557e9e4536cSHong Zhang #if defined(DEBUG_MATMATMULT) 558e9e4536cSHong Zhang if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] C_seq Numeric is done...\n",rank); 5596291f5a6SHong Zhang #endif 5602499ec78SHong Zhang 5612499ec78SHong Zhang /* create mpi matrix C by concatinating C_seq */ 5622499ec78SHong Zhang ierr = PetscObjectReference((PetscObject)mult->C_seq);CHKERRQ(ierr); /* prevent C_seq being destroyed by MatMerge() */ 563*9b8102ccSHong Zhang ierr = MatMergeSymbolic(((PetscObject)A)->comm,mult->C_seq,B->cmap->n,&AB);CHKERRQ(ierr); 564*9b8102ccSHong Zhang ierr = MatMergeNumeric(((PetscObject)A)->comm,mult->C_seq,B->cmap->n,AB);CHKERRQ(ierr); 5656291f5a6SHong Zhang #if defined(DEBUG_MATMATMULT) 5666291f5a6SHong Zhang if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] Merge is done...\n",rank); 5676291f5a6SHong Zhang #endif 5682499ec78SHong Zhang 5692499ec78SHong Zhang /* attach the supporting struct to C for reuse of symbolic C */ 570776b82aeSLisandro Dalcin ierr = PetscContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr); 571776b82aeSLisandro Dalcin ierr = PetscContainerSetPointer(container,mult);CHKERRQ(ierr); 572776b82aeSLisandro Dalcin ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_Mat_MatMatMultMPI);CHKERRQ(ierr); 573d5821860SHong Zhang ierr = PetscObjectCompose((PetscObject)AB,"Mat_MatMatMultMPI",(PetscObject)container);CHKERRQ(ierr); 574b160f555SLisandro Dalcin ierr = PetscContainerDestroy(&container);CHKERRQ(ierr); 5753f7ae0dbSHong Zhang mult->skipNumeric = PETSC_TRUE; /* a numeric product is done here */ 576d5821860SHong Zhang mult->destroy = AB->ops->destroy; 577d5821860SHong Zhang mult->duplicate = AB->ops->duplicate; 578cf1a0672SHong Zhang AB->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ_Scalable_32; 579a1a4e84aSHong Zhang AB->ops->destroy = MatDestroy_MPIAIJ_MatMatMult_32; 580a1a4e84aSHong Zhang AB->ops->duplicate = MatDuplicate_MPIAIJ_MatMatMult_32; 5818cdbd757SHong Zhang AB->ops->matmult = MatMatMult_MPIAIJ_MPIAIJ; 582d5821860SHong Zhang *C = AB; 5832499ec78SHong Zhang PetscFunctionReturn(0); 5842499ec78SHong Zhang } 5852499ec78SHong Zhang 5869123193aSHong Zhang #undef __FUNCT__ 5879123193aSHong Zhang #define __FUNCT__ "MatMatMult_MPIAIJ_MPIDense" 5889123193aSHong Zhang PetscErrorCode MatMatMult_MPIAIJ_MPIDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 5899123193aSHong Zhang { 5909123193aSHong Zhang PetscErrorCode ierr; 5919123193aSHong Zhang 5929123193aSHong Zhang PetscFunctionBegin; 5939123193aSHong Zhang if (scall == MAT_INITIAL_MATRIX){ 5949123193aSHong Zhang ierr = MatMatMultSymbolic_MPIAIJ_MPIDense(A,B,fill,C);CHKERRQ(ierr); 5959123193aSHong Zhang } 5969123193aSHong Zhang ierr = MatMatMultNumeric_MPIAIJ_MPIDense(A,B,*C);CHKERRQ(ierr); 5979123193aSHong Zhang PetscFunctionReturn(0); 5989123193aSHong Zhang } 5999123193aSHong Zhang 6004324174eSBarry Smith typedef struct { 6014324174eSBarry Smith Mat workB; 6024324174eSBarry Smith PetscScalar *rvalues,*svalues; 6034324174eSBarry Smith MPI_Request *rwaits,*swaits; 6044324174eSBarry Smith } MPIAIJ_MPIDense; 6054324174eSBarry Smith 6067af9e4fdSHong Zhang #undef __FUNCT__ 6077af9e4fdSHong Zhang #define __FUNCT__ "MPIAIJ_MPIDenseDestroy" 6084324174eSBarry Smith PetscErrorCode MPIAIJ_MPIDenseDestroy(void *ctx) 6094324174eSBarry Smith { 6104324174eSBarry Smith MPIAIJ_MPIDense *contents = (MPIAIJ_MPIDense*) ctx; 6114324174eSBarry Smith PetscErrorCode ierr; 6124324174eSBarry Smith 6134324174eSBarry Smith PetscFunctionBegin; 6146bf464f9SBarry Smith ierr = MatDestroy(&contents->workB);CHKERRQ(ierr); 6154324174eSBarry Smith ierr = PetscFree4(contents->rvalues,contents->svalues,contents->rwaits,contents->swaits);CHKERRQ(ierr); 6164324174eSBarry Smith ierr = PetscFree(contents);CHKERRQ(ierr); 6174324174eSBarry Smith PetscFunctionReturn(0); 6184324174eSBarry Smith } 6194324174eSBarry Smith 6209123193aSHong Zhang #undef __FUNCT__ 6219123193aSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIDense" 6229123193aSHong Zhang PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIDense(Mat A,Mat B,PetscReal fill,Mat *C) 6239123193aSHong Zhang { 6249123193aSHong Zhang PetscErrorCode ierr; 625f32d5d43SBarry Smith Mat_MPIAIJ *aij = (Mat_MPIAIJ*) A->data; 626d0f46423SBarry Smith PetscInt nz = aij->B->cmap->n; 627bf0cc555SLisandro Dalcin PetscContainer container; 6284324174eSBarry Smith MPIAIJ_MPIDense *contents; 6294324174eSBarry Smith VecScatter ctx = aij->Mvctx; 6304324174eSBarry Smith VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata; 6314324174eSBarry Smith VecScatter_MPI_General *to = ( VecScatter_MPI_General*) ctx->todata; 632d0f46423SBarry Smith PetscInt m=A->rmap->n,n=B->cmap->n; 6339123193aSHong Zhang 6349123193aSHong Zhang PetscFunctionBegin; 635cb2480eaSBarry Smith ierr = MatCreate(((PetscObject)B)->comm,C);CHKERRQ(ierr); 636d0f46423SBarry Smith ierr = MatSetSizes(*C,m,n,A->rmap->N,B->cmap->N);CHKERRQ(ierr); 637cb2480eaSBarry Smith ierr = MatSetType(*C,MATMPIDENSE);CHKERRQ(ierr); 638cb2480eaSBarry Smith ierr = MatAssemblyBegin(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 639cb2480eaSBarry Smith ierr = MatAssemblyEnd(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 6408cdbd757SHong Zhang (*C)->ops->matmult = MatMatMult_MPIAIJ_MPIDense; 641f32d5d43SBarry Smith 6424324174eSBarry Smith ierr = PetscNew(MPIAIJ_MPIDense,&contents);CHKERRQ(ierr); 643f32d5d43SBarry Smith /* Create work matrix used to store off processor rows of B needed for local product */ 644d0f46423SBarry Smith ierr = MatCreateSeqDense(PETSC_COMM_SELF,nz,B->cmap->N,PETSC_NULL,&contents->workB);CHKERRQ(ierr); 6454324174eSBarry Smith /* Create work arrays needed */ 646d0f46423SBarry Smith ierr = PetscMalloc4(B->cmap->N*from->starts[from->n],PetscScalar,&contents->rvalues, 647d0f46423SBarry Smith B->cmap->N*to->starts[to->n],PetscScalar,&contents->svalues, 6484324174eSBarry Smith from->n,MPI_Request,&contents->rwaits, 6494324174eSBarry Smith to->n,MPI_Request,&contents->swaits);CHKERRQ(ierr); 6504324174eSBarry Smith 651bf0cc555SLisandro Dalcin ierr = PetscContainerCreate(((PetscObject)A)->comm,&container);CHKERRQ(ierr); 652bf0cc555SLisandro Dalcin ierr = PetscContainerSetPointer(container,contents);CHKERRQ(ierr); 653bf0cc555SLisandro Dalcin ierr = PetscContainerSetUserDestroy(container,MPIAIJ_MPIDenseDestroy);CHKERRQ(ierr); 654bf0cc555SLisandro Dalcin ierr = PetscObjectCompose((PetscObject)(*C),"workB",(PetscObject)container);CHKERRQ(ierr); 655bf0cc555SLisandro Dalcin ierr = PetscContainerDestroy(&container);CHKERRQ(ierr); 6569123193aSHong Zhang PetscFunctionReturn(0); 6579123193aSHong Zhang } 6589123193aSHong Zhang 6597af9e4fdSHong Zhang #undef __FUNCT__ 6607af9e4fdSHong Zhang #define __FUNCT__ "MatMPIDenseScatter" 661f32d5d43SBarry Smith /* 6622636ff26SBarry Smith Performs an efficient scatter on the rows of B needed by this process; this is 6632636ff26SBarry Smith a modification of the VecScatterBegin_() routines. 664f32d5d43SBarry Smith */ 6654324174eSBarry Smith PetscErrorCode MatMPIDenseScatter(Mat A,Mat B,Mat C,Mat *outworkB) 666f32d5d43SBarry Smith { 667f32d5d43SBarry Smith Mat_MPIAIJ *aij = (Mat_MPIAIJ*)A->data; 668f32d5d43SBarry Smith PetscErrorCode ierr; 669f32d5d43SBarry Smith PetscScalar *b,*w,*svalues,*rvalues; 670f32d5d43SBarry Smith VecScatter ctx = aij->Mvctx; 671f32d5d43SBarry Smith VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata; 672f32d5d43SBarry Smith VecScatter_MPI_General *to = ( VecScatter_MPI_General*) ctx->todata; 673f32d5d43SBarry Smith PetscInt i,j,k; 674aa5bb8c0SSatish Balay PetscInt *sindices,*sstarts,*rindices,*rstarts; 675aa5bb8c0SSatish Balay PetscMPIInt *sprocs,*rprocs,nrecvs; 676f32d5d43SBarry Smith MPI_Request *swaits,*rwaits; 6777adad957SLisandro Dalcin MPI_Comm comm = ((PetscObject)A)->comm; 678d0f46423SBarry Smith PetscMPIInt tag = ((PetscObject)ctx)->tag,ncols = B->cmap->N, nrows = aij->B->cmap->n,imdex,nrowsB = B->rmap->n; 679f32d5d43SBarry Smith MPI_Status status; 6804324174eSBarry Smith MPIAIJ_MPIDense *contents; 681bf0cc555SLisandro Dalcin PetscContainer container; 6824324174eSBarry Smith Mat workB; 683f32d5d43SBarry Smith 684f32d5d43SBarry Smith PetscFunctionBegin; 685bf0cc555SLisandro Dalcin ierr = PetscObjectQuery((PetscObject)C,"workB",(PetscObject*)&container);CHKERRQ(ierr); 68629825010SHong Zhang if (!container) SETERRQ(comm,PETSC_ERR_PLIB,"Container does not exist"); 687bf0cc555SLisandro Dalcin ierr = PetscContainerGetPointer(container,(void**)&contents);CHKERRQ(ierr); 6884324174eSBarry Smith 6894324174eSBarry Smith workB = *outworkB = contents->workB; 690cf1a0672SHong 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); 691f32d5d43SBarry Smith sindices = to->indices; 692f32d5d43SBarry Smith sstarts = to->starts; 693f32d5d43SBarry Smith sprocs = to->procs; 6944324174eSBarry Smith swaits = contents->swaits; 6954324174eSBarry Smith svalues = contents->svalues; 696f32d5d43SBarry Smith 697f32d5d43SBarry Smith rindices = from->indices; 698f32d5d43SBarry Smith rstarts = from->starts; 699f32d5d43SBarry Smith rprocs = from->procs; 7004324174eSBarry Smith rwaits = contents->rwaits; 7014324174eSBarry Smith rvalues = contents->rvalues; 702f32d5d43SBarry Smith 703f32d5d43SBarry Smith ierr = MatGetArray(B,&b);CHKERRQ(ierr); 704f32d5d43SBarry Smith ierr = MatGetArray(workB,&w);CHKERRQ(ierr); 705f32d5d43SBarry Smith 706f32d5d43SBarry Smith for (i=0; i<from->n; i++) { 707f32d5d43SBarry Smith ierr = MPI_Irecv(rvalues+ncols*rstarts[i],ncols*(rstarts[i+1]-rstarts[i]),MPIU_SCALAR,rprocs[i],tag,comm,rwaits+i);CHKERRQ(ierr); 708f32d5d43SBarry Smith } 709f32d5d43SBarry Smith 710f32d5d43SBarry Smith for (i=0; i<to->n; i++) { 711f32d5d43SBarry Smith /* pack a message at a time */ 712f32d5d43SBarry Smith CHKMEMQ; 713f32d5d43SBarry Smith for (j=0; j<sstarts[i+1]-sstarts[i]; j++){ 714f32d5d43SBarry Smith for (k=0; k<ncols; k++) { 715f32d5d43SBarry Smith svalues[ncols*(sstarts[i] + j) + k] = b[sindices[sstarts[i]+j] + nrowsB*k]; 7162636ff26SBarry Smith } 7172636ff26SBarry Smith } 718f32d5d43SBarry Smith CHKMEMQ; 719f32d5d43SBarry Smith ierr = MPI_Isend(svalues+ncols*sstarts[i],ncols*(sstarts[i+1]-sstarts[i]),MPIU_SCALAR,sprocs[i],tag,comm,swaits+i);CHKERRQ(ierr); 720f32d5d43SBarry Smith } 7212636ff26SBarry Smith 722f32d5d43SBarry Smith nrecvs = from->n; 723f32d5d43SBarry Smith while (nrecvs) { 724f32d5d43SBarry Smith ierr = MPI_Waitany(from->n,rwaits,&imdex,&status);CHKERRQ(ierr); 725f32d5d43SBarry Smith nrecvs--; 726f32d5d43SBarry Smith /* unpack a message at a time */ 727f32d5d43SBarry Smith CHKMEMQ; 728f32d5d43SBarry Smith for (j=0; j<rstarts[imdex+1]-rstarts[imdex]; j++){ 729f32d5d43SBarry Smith for (k=0; k<ncols; k++) { 730f32d5d43SBarry Smith w[rindices[rstarts[imdex]+j] + nrows*k] = rvalues[ncols*(rstarts[imdex] + j) + k]; 7312636ff26SBarry Smith } 7322636ff26SBarry Smith } 733f32d5d43SBarry Smith CHKMEMQ; 734f32d5d43SBarry Smith } 735cb9801acSJed Brown if (to->n) {ierr = MPI_Waitall(to->n,swaits,to->sstatus);CHKERRQ(ierr);} 736f32d5d43SBarry Smith 737f32d5d43SBarry Smith ierr = MatRestoreArray(B,&b);CHKERRQ(ierr); 738f32d5d43SBarry Smith ierr = MatRestoreArray(workB,&w);CHKERRQ(ierr); 739f32d5d43SBarry Smith ierr = MatAssemblyBegin(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 740f32d5d43SBarry Smith ierr = MatAssemblyEnd(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 741f32d5d43SBarry Smith PetscFunctionReturn(0); 742f32d5d43SBarry Smith } 743f32d5d43SBarry Smith extern PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat,Mat,Mat); 744f32d5d43SBarry Smith 7459123193aSHong Zhang #undef __FUNCT__ 7469123193aSHong Zhang #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIDense" 7479123193aSHong Zhang PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIDense(Mat A,Mat B,Mat C) 7489123193aSHong Zhang { 7499123193aSHong Zhang PetscErrorCode ierr; 750f32d5d43SBarry Smith Mat_MPIAIJ *aij = (Mat_MPIAIJ*)A->data; 751f32d5d43SBarry Smith Mat_MPIDense *bdense = (Mat_MPIDense*)B->data; 752f32d5d43SBarry Smith Mat_MPIDense *cdense = (Mat_MPIDense*)C->data; 753f32d5d43SBarry Smith Mat workB; 7549123193aSHong Zhang 7559123193aSHong Zhang PetscFunctionBegin; 7569123193aSHong Zhang 757f32d5d43SBarry Smith /* diagonal block of A times all local rows of B*/ 758f32d5d43SBarry Smith ierr = MatMatMultNumeric_SeqAIJ_SeqDense(aij->A,bdense->A,cdense->A);CHKERRQ(ierr); 759f32d5d43SBarry Smith 760f32d5d43SBarry Smith /* get off processor parts of B needed to complete the product */ 7614324174eSBarry Smith ierr = MatMPIDenseScatter(A,B,C,&workB);CHKERRQ(ierr); 762f32d5d43SBarry Smith 763f32d5d43SBarry Smith /* off-diagonal block of A times nonlocal rows of B */ 764f32d5d43SBarry Smith ierr = MatMatMultNumericAdd_SeqAIJ_SeqDense(aij->B,workB,cdense->A);CHKERRQ(ierr); 7659123193aSHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 7669123193aSHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 7679123193aSHong Zhang PetscFunctionReturn(0); 7689123193aSHong Zhang } 769cf1a0672SHong Zhang 7701634b032SHong Zhang #undef __FUNCT__ 771cf1a0672SHong Zhang #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIAIJ_Scalable" 772cf1a0672SHong Zhang PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ_Scalable(Mat A,Mat P,Mat C) 7731634b032SHong Zhang { 7741634b032SHong Zhang PetscErrorCode ierr; 775cf1a0672SHong Zhang Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data,*c=(Mat_MPIAIJ*)C->data; 776cf1a0672SHong Zhang Mat_SeqAIJ *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data; 777cf1a0672SHong Zhang Mat_SeqAIJ *cd=(Mat_SeqAIJ*)(c->A)->data,*co=(Mat_SeqAIJ*)(c->B)->data; 778cf1a0672SHong Zhang PetscInt *adi=ad->i,*adj,*aoi=ao->i,*aoj; 779cf1a0672SHong Zhang PetscScalar *ada,*aoa,*cda=cd->a,*coa=co->a; 780cf1a0672SHong Zhang Mat_SeqAIJ *p_loc,*p_oth; 781cf1a0672SHong Zhang PetscInt *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pj; 78229825010SHong Zhang PetscScalar *pa_loc,*pa_oth,*pa,valtmp,*ca; 78329825010SHong Zhang PetscInt cm=C->rmap->n,anz,pnz; 784cf1a0672SHong Zhang Mat_PtAPMPI *ptap=c->ptap; 78529825010SHong Zhang PetscScalar *apa_sparse=ptap->apa; 786cf1a0672SHong Zhang PetscInt *api,*apj,*apJ,i,j,k,row; 78729825010SHong Zhang PetscInt cstart=C->cmap->rstart; 78829825010SHong Zhang PetscInt cdnz,conz,k0,k1,nextp; 7891634b032SHong Zhang #if defined(DEBUG_MATMATMULT) 790cf1a0672SHong Zhang PetscMPIInt rank=a->rank; 7911634b032SHong Zhang #endif 7921634b032SHong Zhang 7931634b032SHong Zhang PetscFunctionBegin; 7941634b032SHong Zhang #if defined(DEBUG_MATMATMULT) 795cf1a0672SHong Zhang if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] call MatMatMultNumeric_MPIAIJ_MPIAIJ_Scalable()...\n",rank); 7961634b032SHong Zhang #endif 797cf1a0672SHong Zhang 798cf1a0672SHong Zhang /* 1) get P_oth = ptap->P_oth and P_loc = ptap->P_loc */ 799cf1a0672SHong Zhang /*-----------------------------------------------------*/ 800cf1a0672SHong Zhang /* update numerical values of P_oth and P_loc */ 801b7f45c76SHong Zhang ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 802cf1a0672SHong Zhang ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 803cf1a0672SHong Zhang #if defined(DEBUG_MATMATMULT) 804cf1a0672SHong Zhang if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] got P_oth and P_loc...\n",rank); 805cf1a0672SHong Zhang #endif 806cf1a0672SHong Zhang 807cf1a0672SHong Zhang /* 2) compute numeric C_loc = A_loc*P = Ad*P_loc + Ao*P_oth */ 808cf1a0672SHong Zhang /*----------------------------------------------------------*/ 809cf1a0672SHong Zhang /* get data from symbolic products */ 810cf1a0672SHong Zhang p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 811cf1a0672SHong Zhang p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 812cf1a0672SHong Zhang pi_loc=p_loc->i; pj_loc=p_loc->j; pa_loc=p_loc->a; 813cf1a0672SHong Zhang pi_oth=p_oth->i; pj_oth=p_oth->j; pa_oth=p_oth->a; 814cf1a0672SHong Zhang 81529825010SHong Zhang api = ptap->api; 81629825010SHong Zhang apj = ptap->apj; 817cf1a0672SHong Zhang for (i=0; i<cm; i++) { 81829825010SHong Zhang apJ = apj + api[i]; 81929825010SHong Zhang 820cf1a0672SHong Zhang /* diagonal portion of A */ 821cf1a0672SHong Zhang anz = adi[i+1] - adi[i]; 822cf1a0672SHong Zhang adj = ad->j + adi[i]; 823cf1a0672SHong Zhang ada = ad->a + adi[i]; 824cf1a0672SHong Zhang for (j=0; j<anz; j++) { 825cf1a0672SHong Zhang row = adj[j]; 826cf1a0672SHong Zhang pnz = pi_loc[row+1] - pi_loc[row]; 827cf1a0672SHong Zhang pj = pj_loc + pi_loc[row]; 828cf1a0672SHong Zhang pa = pa_loc + pi_loc[row]; 82929825010SHong Zhang /* perform sparse axpy */ 830cf1a0672SHong Zhang valtmp = ada[j]; 83129825010SHong Zhang nextp = 0; 83229825010SHong Zhang for (k=0; nextp<pnz; k++) { 83329825010SHong Zhang if (apJ[k] == pj[nextp]) { /* column of AP == column of P */ 83429825010SHong Zhang apa_sparse[k] += valtmp*pa[nextp++]; 83529825010SHong Zhang } 8361634b032SHong Zhang } 837cf1a0672SHong Zhang ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr); 838cf1a0672SHong Zhang } 8391634b032SHong Zhang 840cf1a0672SHong Zhang /* off-diagonal portion of A */ 841cf1a0672SHong Zhang anz = aoi[i+1] - aoi[i]; 842cf1a0672SHong Zhang aoj = ao->j + aoi[i]; 843cf1a0672SHong Zhang aoa = ao->a + aoi[i]; 844cf1a0672SHong Zhang for (j=0; j<anz; j++) { 845cf1a0672SHong Zhang row = aoj[j]; 846cf1a0672SHong Zhang pnz = pi_oth[row+1] - pi_oth[row]; 847cf1a0672SHong Zhang pj = pj_oth + pi_oth[row]; 848cf1a0672SHong Zhang pa = pa_oth + pi_oth[row]; 84929825010SHong Zhang /* perform sparse axpy */ 850cf1a0672SHong Zhang valtmp = aoa[j]; 85129825010SHong Zhang nextp = 0; 85229825010SHong Zhang for (k=0; nextp<pnz; k++) { 85329825010SHong Zhang if (apJ[k] == pj[nextp]) { /* column of AP == column of P */ 85429825010SHong Zhang apa_sparse[k] += valtmp*pa[nextp++]; 85529825010SHong Zhang } 856cf1a0672SHong Zhang } 857cf1a0672SHong Zhang ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr); 858cf1a0672SHong Zhang } 8591634b032SHong Zhang 860cf1a0672SHong Zhang /* set values in C */ 861cf1a0672SHong Zhang cdnz = cd->i[i+1] - cd->i[i]; 862cf1a0672SHong Zhang conz = co->i[i+1] - co->i[i]; 8631634b032SHong Zhang 864cf1a0672SHong Zhang /* 1st off-diagoanl part of C */ 865cf1a0672SHong Zhang ca = coa + co->i[i]; 866cf1a0672SHong Zhang k = 0; 867cf1a0672SHong Zhang for (k0=0; k0<conz; k0++){ 868cf1a0672SHong Zhang if (apJ[k] >= cstart) break; 86929825010SHong Zhang ca[k0] = apa_sparse[k]; 87029825010SHong Zhang apa_sparse[k] = 0.0; 871cf1a0672SHong Zhang k++; 872cf1a0672SHong Zhang } 8731634b032SHong Zhang 874cf1a0672SHong Zhang /* diagonal part of C */ 875cf1a0672SHong Zhang ca = cda + cd->i[i]; 876cf1a0672SHong Zhang for (k1=0; k1<cdnz; k1++){ 87729825010SHong Zhang ca[k1] = apa_sparse[k]; 87829825010SHong Zhang apa_sparse[k] = 0.0; 879cf1a0672SHong Zhang k++; 880cf1a0672SHong Zhang } 881cf1a0672SHong Zhang 882cf1a0672SHong Zhang /* 2nd off-diagoanl part of C */ 883cf1a0672SHong Zhang ca = coa + co->i[i]; 884cf1a0672SHong Zhang for (; k0<conz; k0++){ 88529825010SHong Zhang ca[k0] = apa_sparse[k]; 88629825010SHong Zhang apa_sparse[k] = 0.0; 887cf1a0672SHong Zhang k++; 888cf1a0672SHong Zhang } 889cf1a0672SHong Zhang } 890cf1a0672SHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 891cf1a0672SHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 892cf1a0672SHong Zhang PetscFunctionReturn(0); 893cf1a0672SHong Zhang } 894cf1a0672SHong Zhang 895cf1a0672SHong Zhang /* same as MatMatMultSymbolic_MPIAIJ_MPIAIJ(), except using LLCondensed to avoid O(BN) memory requirement */ 896cf1a0672SHong Zhang #undef __FUNCT__ 897cf1a0672SHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable" 898cf1a0672SHong Zhang PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable(Mat A,Mat P,PetscReal fill,Mat *C) 899cf1a0672SHong Zhang { 900cf1a0672SHong Zhang PetscErrorCode ierr; 901cf1a0672SHong Zhang MPI_Comm comm=((PetscObject)A)->comm; 902cf1a0672SHong Zhang Mat Cmpi; 903cf1a0672SHong Zhang Mat_PtAPMPI *ptap; 904cf1a0672SHong Zhang PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 905cf1a0672SHong Zhang Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data,*c; 906cf1a0672SHong Zhang Mat_SeqAIJ *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_loc,*p_oth; 907cf1a0672SHong Zhang PetscInt *pi_loc,*pj_loc,*pi_oth,*pj_oth,*dnz,*onz; 908cf1a0672SHong Zhang PetscInt *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,rstart=A->rmap->rstart; 909f84c536bSHong Zhang PetscInt i,pnz,row,*api,*apj,*Jptr,apnz,nspacedouble=0,j,nzi,*lnk,apnz_max=0; 910cf1a0672SHong Zhang PetscInt am=A->rmap->n,pN=P->cmap->N,pn=P->cmap->n,pm=P->rmap->n; 911cf1a0672SHong Zhang PetscInt nlnk_max,armax,prmax; 912cf1a0672SHong Zhang PetscReal afill; 913cf1a0672SHong Zhang PetscScalar *apa; 914cf1a0672SHong Zhang #if defined(DEBUG_MATMATMULT) 915cf1a0672SHong Zhang PetscMPIInt rank=a->rank; 916cf1a0672SHong Zhang #endif 917cf1a0672SHong Zhang 918cf1a0672SHong Zhang PetscFunctionBegin; 919cf1a0672SHong Zhang #if defined(DEBUG_MATMATMULT) 920cf1a0672SHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 921cf1a0672SHong Zhang if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] call MatMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable()...\n",rank); 922cf1a0672SHong Zhang #endif 923cf1a0672SHong Zhang 924cf1a0672SHong Zhang /* create struct Mat_PtAPMPI and attached it to C later */ 925cf1a0672SHong Zhang ierr = PetscNew(Mat_PtAPMPI,&ptap);CHKERRQ(ierr); 926cf1a0672SHong Zhang 927cf1a0672SHong Zhang /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */ 928b7f45c76SHong Zhang ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 929cf1a0672SHong Zhang #if defined(DEBUG_MATMATMULT) 930cf1a0672SHong Zhang if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] P_oth is done...\n",rank); 931cf1a0672SHong Zhang #endif 932cf1a0672SHong Zhang /* get P_loc by taking all local rows of P */ 933cf1a0672SHong Zhang ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 934cf1a0672SHong Zhang 935cf1a0672SHong Zhang p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 936cf1a0672SHong Zhang p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 937cf1a0672SHong Zhang pi_loc = p_loc->i; pj_loc = p_loc->j; 938cf1a0672SHong Zhang pi_oth = p_oth->i; pj_oth = p_oth->j; 939cf1a0672SHong Zhang 940cf1a0672SHong Zhang #if defined(DEBUG_MATMATMULT) 941cf1a0672SHong 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); 942cf1a0672SHong Zhang #endif 943cf1a0672SHong Zhang 944cf1a0672SHong Zhang /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */ 945cf1a0672SHong Zhang /*-------------------------------------------------------------------*/ 946cf1a0672SHong Zhang ierr = PetscMalloc((am+2)*sizeof(PetscInt),&api);CHKERRQ(ierr); 947cf1a0672SHong Zhang ptap->api = api; 948cf1a0672SHong Zhang api[0] = 0; 949cf1a0672SHong Zhang 950cf1a0672SHong Zhang /* create and initialize a linked list */ 951cf1a0672SHong Zhang armax = ad->rmax+ao->rmax; 952cf1a0672SHong Zhang prmax = PetscMax(p_loc->rmax,p_oth->rmax); 953cf1a0672SHong Zhang nlnk_max = armax*prmax; 954cf1a0672SHong Zhang if (!nlnk_max || nlnk_max > pN) nlnk_max = pN; 955cf1a0672SHong Zhang #if defined(DEBUG_MATMATMULT) 956cf1a0672SHong 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); 957cf1a0672SHong Zhang #endif 958f84c536bSHong Zhang ierr = PetscLLCondensedCreate_Scalable(nlnk_max,&lnk);CHKERRQ(ierr); 959cf1a0672SHong Zhang 960cf1a0672SHong Zhang /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */ 961cf1a0672SHong Zhang ierr = PetscFreeSpaceGet((PetscInt)(fill*(adi[am]+aoi[am]+pi_loc[pm])),&free_space);CHKERRQ(ierr); 962cf1a0672SHong Zhang current_space = free_space; 963cf1a0672SHong Zhang 964cf1a0672SHong Zhang ierr = MatPreallocateInitialize(comm,am,pn,dnz,onz);CHKERRQ(ierr); 965cf1a0672SHong Zhang for (i=0; i<am; i++) { 966cf1a0672SHong Zhang apnz = 0; 967cf1a0672SHong Zhang /* diagonal portion of A */ 968cf1a0672SHong Zhang nzi = adi[i+1] - adi[i]; 969cf1a0672SHong Zhang for (j=0; j<nzi; j++){ 970cf1a0672SHong Zhang row = *adj++; 971cf1a0672SHong Zhang pnz = pi_loc[row+1] - pi_loc[row]; 972cf1a0672SHong Zhang Jptr = pj_loc + pi_loc[row]; 973cf1a0672SHong Zhang /* add non-zero cols of P into the sorted linked list lnk */ 974f84c536bSHong Zhang ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr); 975cf1a0672SHong Zhang } 976cf1a0672SHong Zhang /* off-diagonal portion of A */ 977cf1a0672SHong Zhang nzi = aoi[i+1] - aoi[i]; 978cf1a0672SHong Zhang for (j=0; j<nzi; j++){ 979cf1a0672SHong Zhang row = *aoj++; 980cf1a0672SHong Zhang pnz = pi_oth[row+1] - pi_oth[row]; 981cf1a0672SHong Zhang Jptr = pj_oth + pi_oth[row]; 982f84c536bSHong Zhang ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr); 983cf1a0672SHong Zhang } 984cf1a0672SHong Zhang 985f84c536bSHong Zhang apnz = *lnk; 986cf1a0672SHong Zhang api[i+1] = api[i] + apnz; 987cf1a0672SHong Zhang if (apnz > apnz_max) apnz_max = apnz; 988cf1a0672SHong Zhang 989cf1a0672SHong Zhang /* if free space is not available, double the total space in the list */ 990cf1a0672SHong Zhang if (current_space->local_remaining<apnz) { 991cf1a0672SHong Zhang ierr = PetscFreeSpaceGet(apnz+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 992cf1a0672SHong Zhang nspacedouble++; 993cf1a0672SHong Zhang } 994cf1a0672SHong Zhang 995cf1a0672SHong Zhang /* Copy data into free space, then initialize lnk */ 996f84c536bSHong Zhang ierr = PetscLLCondensedClean_Scalable(apnz,current_space->array,lnk);CHKERRQ(ierr); 997cf1a0672SHong Zhang ierr = MatPreallocateSet(i+rstart,apnz,current_space->array,dnz,onz);CHKERRQ(ierr); 998cf1a0672SHong Zhang current_space->array += apnz; 999cf1a0672SHong Zhang current_space->local_used += apnz; 1000cf1a0672SHong Zhang current_space->local_remaining -= apnz; 1001cf1a0672SHong Zhang } 1002cf1a0672SHong Zhang 1003cf1a0672SHong Zhang /* Allocate space for apj, initialize apj, and */ 1004cf1a0672SHong Zhang /* destroy list of free space and other temporary array(s) */ 1005cf1a0672SHong Zhang ierr = PetscMalloc((api[am]+1)*sizeof(PetscInt),&ptap->apj);CHKERRQ(ierr); 1006cf1a0672SHong Zhang apj = ptap->apj; 1007cf1a0672SHong Zhang ierr = PetscFreeSpaceContiguous(&free_space,ptap->apj);CHKERRQ(ierr); 1008f84c536bSHong Zhang ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr); 1009cf1a0672SHong Zhang #if defined(DEBUG_MATMATMULT) 1010cf1a0672SHong Zhang if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] AP is done..., apnz_max %d\n",rank,apnz_max); 1011cf1a0672SHong Zhang #endif 1012cf1a0672SHong Zhang 1013cf1a0672SHong Zhang /* create and assemble symbolic parallel matrix Cmpi */ 1014cf1a0672SHong Zhang /*----------------------------------------------------*/ 1015cf1a0672SHong Zhang ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr); 1016cf1a0672SHong Zhang ierr = MatSetSizes(Cmpi,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 1017cf1a0672SHong Zhang ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr); 1018cf1a0672SHong Zhang ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr); 1019cf1a0672SHong Zhang ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 1020cf1a0672SHong Zhang ierr = MatSetBlockSize(Cmpi,1);CHKERRQ(ierr); 1021cf1a0672SHong Zhang 102229825010SHong Zhang /* malloc apa for assembly Cmpi */ 1023cf1a0672SHong Zhang ierr = PetscMalloc(apnz_max*sizeof(PetscScalar),&apa);CHKERRQ(ierr); 1024cf1a0672SHong Zhang ierr = PetscMemzero(apa,apnz_max*sizeof(PetscScalar));CHKERRQ(ierr); 102529825010SHong Zhang ptap->apa = apa; 1026cf1a0672SHong Zhang for (i=0; i<am; i++){ 1027cf1a0672SHong Zhang row = i + rstart; 1028cf1a0672SHong Zhang apnz = api[i+1] - api[i]; 1029cf1a0672SHong Zhang ierr = MatSetValues(Cmpi,1,&row,apnz,apj,apa,INSERT_VALUES);CHKERRQ(ierr); 1030cf1a0672SHong Zhang apj += apnz; 1031cf1a0672SHong Zhang } 1032cf1a0672SHong Zhang ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1033cf1a0672SHong Zhang ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1034cf1a0672SHong Zhang 1035cf1a0672SHong Zhang ptap->destroy = Cmpi->ops->destroy; 1036cf1a0672SHong Zhang ptap->duplicate = Cmpi->ops->duplicate; 1037cf1a0672SHong Zhang Cmpi->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ_Scalable; 1038cf1a0672SHong Zhang Cmpi->ops->destroy = MatDestroy_MPIAIJ_MatMatMult; 1039cf1a0672SHong Zhang Cmpi->ops->duplicate = MatDuplicate_MPIAIJ_MatMatMult; 1040cf1a0672SHong Zhang 1041cf1a0672SHong Zhang /* attach the supporting struct to Cmpi for reuse */ 1042cf1a0672SHong Zhang c = (Mat_MPIAIJ*)Cmpi->data; 1043cf1a0672SHong Zhang c->ptap = ptap; 1044cf1a0672SHong Zhang 1045cf1a0672SHong Zhang *C = Cmpi; 1046cf1a0672SHong Zhang 1047cf1a0672SHong Zhang /* set MatInfo */ 1048cf1a0672SHong Zhang afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]) + 1.e-5; 1049cf1a0672SHong Zhang if (afill < 1.0) afill = 1.0; 1050cf1a0672SHong Zhang Cmpi->info.mallocs = nspacedouble; 1051cf1a0672SHong Zhang Cmpi->info.fill_ratio_given = fill; 1052cf1a0672SHong Zhang Cmpi->info.fill_ratio_needed = afill; 1053cf1a0672SHong Zhang 1054cf1a0672SHong Zhang #if defined(PETSC_USE_INFO) 1055cf1a0672SHong Zhang if (api[am]) { 1056cf1a0672SHong Zhang ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr); 1057cf1a0672SHong Zhang ierr = PetscInfo1(Cmpi,"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr); 1058cf1a0672SHong Zhang } else { 1059cf1a0672SHong Zhang ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr); 1060cf1a0672SHong Zhang } 1061cf1a0672SHong Zhang #endif 10621634b032SHong Zhang PetscFunctionReturn(0); 10631634b032SHong Zhang } 1064f96d37fcSHong Zhang 1065f96d37fcSHong Zhang /*-------------------------------------------------------------------------*/ 1066f96d37fcSHong Zhang #undef __FUNCT__ 1067f96d37fcSHong Zhang #define __FUNCT__ "MatTransposeMatMult_MPIAIJ_MPIAIJ" 1068c5af039cSHong Zhang PetscErrorCode MatTransposeMatMult_MPIAIJ_MPIAIJ(Mat P,Mat A,MatReuse scall,PetscReal fill,Mat *C) 1069f96d37fcSHong Zhang { 1070f96d37fcSHong Zhang PetscErrorCode ierr; 1071d6ab1dc0SHong Zhang PetscBool scalable=PETSC_FALSE; 1072f96d37fcSHong Zhang 1073f96d37fcSHong Zhang PetscFunctionBegin; 1074f96d37fcSHong Zhang if (scall == MAT_INITIAL_MATRIX){ 1075d6ab1dc0SHong Zhang ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr); 1076d6ab1dc0SHong Zhang ierr = PetscOptionsBool("-mattransposematmult_scalable","Use a scalable but slower C=Pt*A","",scalable,&scalable,PETSC_NULL);CHKERRQ(ierr); 1077d6ab1dc0SHong Zhang if (scalable){ 1078d6ab1dc0SHong Zhang ierr = MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable(P,A,fill,C);CHKERRQ(ierr); 1079d6ab1dc0SHong Zhang } else { 1080c5af039cSHong Zhang ierr = MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(P,A,fill,C);CHKERRQ(ierr); 1081f96d37fcSHong Zhang } 1082d6ab1dc0SHong Zhang ierr = PetscOptionsEnd();CHKERRQ(ierr); 1083d6ab1dc0SHong Zhang } 1084d6ab1dc0SHong Zhang ierr = (*(*C)->ops->mattransposemultnumeric)(P,A,*C);CHKERRQ(ierr); 1085f96d37fcSHong Zhang PetscFunctionReturn(0); 1086f96d37fcSHong Zhang } 1087f96d37fcSHong Zhang 1088f96d37fcSHong Zhang #undef __FUNCT__ 1089f96d37fcSHong Zhang #define __FUNCT__ "MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ" 1090c5af039cSHong Zhang PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ(Mat P,Mat A,Mat C) 1091f96d37fcSHong Zhang { 1092c5af039cSHong Zhang PetscErrorCode ierr; 1093c5af039cSHong Zhang Mat_Merge_SeqsToMPI *merge; 1094497f5370SHong Zhang Mat_MPIAIJ *p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data; 1095c5af039cSHong Zhang Mat_SeqAIJ *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data; 1096c5af039cSHong Zhang Mat_PtAPMPI *ptap; 1097d6ab1dc0SHong Zhang PetscInt *adj,*aJ; 1098497f5370SHong Zhang PetscInt i,j,k,anz,pnz,row,*cj; 1099e745005bSHong Zhang MatScalar *ada,*aval,*ca,valtmp; 1100c5af039cSHong Zhang PetscInt am=A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n; 1101c5af039cSHong Zhang MPI_Comm comm=((PetscObject)C)->comm; 1102c5af039cSHong Zhang PetscMPIInt size,rank,taga,*len_s; 1103c5af039cSHong Zhang PetscInt *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci; 1104c5af039cSHong Zhang PetscInt **buf_ri,**buf_rj; 1105c5af039cSHong Zhang PetscInt cnz=0,*bj_i,*bi,*bj,bnz,nextcj; /* bi,bj,ba: local array of C(mpi mat) */ 1106c5af039cSHong Zhang MPI_Request *s_waits,*r_waits; 1107c5af039cSHong Zhang MPI_Status *status; 1108c5af039cSHong Zhang MatScalar **abuf_r,*ba_i,*pA,*coa,*ba; 1109d6ab1dc0SHong Zhang PetscInt *ai,*aj,*coi,*coj; 1110497f5370SHong Zhang PetscInt *poJ=po->j,*pdJ=pd->j; 1111c5af039cSHong Zhang Mat A_loc; 1112c5af039cSHong Zhang Mat_SeqAIJ *a_loc; 1113f96d37fcSHong Zhang 1114f96d37fcSHong Zhang PetscFunctionBegin; 1115c5af039cSHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1116c5af039cSHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 1117c5af039cSHong Zhang 1118c5af039cSHong Zhang ptap = c->ptap; 1119c5af039cSHong Zhang merge = ptap->merge; 1120c5af039cSHong Zhang 1121d6ab1dc0SHong Zhang #if defined(DEBUG_MATTrMatMult) 1122d6ab1dc0SHong Zhang ierr = PetscSynchronizedPrintf(comm,"[%d] TransposeMatMultNumeric: Crmax %d \n",rank,ptap->rmax); 1123d6ab1dc0SHong Zhang ierr = PetscSynchronizedFlush(comm);CHKERRQ(ierr); 1124d6ab1dc0SHong Zhang #endif 1125d6ab1dc0SHong Zhang 1126c5af039cSHong Zhang /* 2) compute numeric C_seq = P_loc^T*A_loc*P - dominating part */ 1127c5af039cSHong Zhang /*--------------------------------------------------------------*/ 1128c5af039cSHong Zhang /* get data from symbolic products */ 1129c5af039cSHong Zhang coi = merge->coi; coj = merge->coj; 1130c5af039cSHong Zhang ierr = PetscMalloc((coi[pon]+1)*sizeof(MatScalar),&coa);CHKERRQ(ierr); 1131c5af039cSHong Zhang ierr = PetscMemzero(coa,coi[pon]*sizeof(MatScalar));CHKERRQ(ierr); 1132c5af039cSHong Zhang 1133c5af039cSHong Zhang bi = merge->bi; bj = merge->bj; 1134c5af039cSHong Zhang owners = merge->rowmap->range; 1135c5af039cSHong Zhang ierr = PetscMalloc((bi[cm]+1)*sizeof(MatScalar),&ba);CHKERRQ(ierr); 1136c5af039cSHong Zhang ierr = PetscMemzero(ba,bi[cm]*sizeof(MatScalar));CHKERRQ(ierr); 1137c5af039cSHong Zhang 1138c5af039cSHong Zhang /* get A_loc by taking all local rows of A */ 1139c5af039cSHong Zhang A_loc = ptap->A_loc; 1140c5af039cSHong Zhang ierr = MatMPIAIJGetLocalMat(A,MAT_REUSE_MATRIX,&A_loc);CHKERRQ(ierr); 1141c5af039cSHong Zhang a_loc = (Mat_SeqAIJ*)(A_loc)->data; 1142d6ab1dc0SHong Zhang ai = a_loc->i; 1143d6ab1dc0SHong Zhang aj = a_loc->j; 1144c5af039cSHong Zhang 1145e745005bSHong Zhang ierr = PetscMalloc((A->cmap->N)*sizeof(PetscScalar),&aval);CHKERRQ(ierr); /* non-scalable!!! */ 1146e745005bSHong Zhang ierr = PetscMemzero(aval,A->cmap->N*sizeof(PetscScalar));CHKERRQ(ierr); 1147c5af039cSHong Zhang 1148c5af039cSHong Zhang for (i=0; i<am; i++) { 1149e745005bSHong Zhang /* 2-a) put A[i,:] to dense array aval */ 1150d6ab1dc0SHong Zhang anz = ai[i+1] - ai[i]; 1151d6ab1dc0SHong Zhang adj = aj + ai[i]; 1152d6ab1dc0SHong Zhang ada = a_loc->a + ai[i]; 1153c5af039cSHong Zhang for (j=0; j<anz; j++){ 1154e745005bSHong Zhang aval[adj[j]] = ada[j]; 1155c5af039cSHong Zhang } 1156c5af039cSHong Zhang 1157c5af039cSHong Zhang /* 2-b) Compute Cseq = P_loc[i,:]^T*A[i,:] using outer product */ 1158c5af039cSHong Zhang /*--------------------------------------------------------------*/ 1159c5af039cSHong Zhang /* put the value into Co=(p->B)^T*A (off-diagonal part, send to others) */ 1160c5af039cSHong Zhang pnz = po->i[i+1] - po->i[i]; 1161c5af039cSHong Zhang poJ = po->j + po->i[i]; 1162c5af039cSHong Zhang pA = po->a + po->i[i]; 1163c5af039cSHong Zhang for (j=0; j<pnz; j++){ 1164c5af039cSHong Zhang row = poJ[j]; 1165c5af039cSHong Zhang cnz = coi[row+1] - coi[row]; 1166c5af039cSHong Zhang cj = coj + coi[row]; 1167c5af039cSHong Zhang ca = coa + coi[row]; 1168c5af039cSHong Zhang /* perform dense axpy */ 1169c5af039cSHong Zhang valtmp = pA[j]; 1170c5af039cSHong Zhang for (k=0; k<cnz; k++) { 1171e745005bSHong Zhang ca[k] += valtmp*aval[cj[k]]; 1172c5af039cSHong Zhang } 1173c5af039cSHong Zhang ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr); 1174c5af039cSHong Zhang } 1175c5af039cSHong Zhang 1176c5af039cSHong Zhang /* put the value into Cd (diagonal part) */ 1177c5af039cSHong Zhang pnz = pd->i[i+1] - pd->i[i]; 1178c5af039cSHong Zhang pdJ = pd->j + pd->i[i]; 1179c5af039cSHong Zhang pA = pd->a + pd->i[i]; 1180c5af039cSHong Zhang for (j=0; j<pnz; j++){ 1181c5af039cSHong Zhang row = pdJ[j]; 1182c5af039cSHong Zhang cnz = bi[row+1] - bi[row]; 1183c5af039cSHong Zhang cj = bj + bi[row]; 1184c5af039cSHong Zhang ca = ba + bi[row]; 1185c5af039cSHong Zhang /* perform dense axpy */ 1186c5af039cSHong Zhang valtmp = pA[j]; 1187c5af039cSHong Zhang for (k=0; k<cnz; k++) { 1188e745005bSHong Zhang ca[k] += valtmp*aval[cj[k]]; 1189c5af039cSHong Zhang } 1190c5af039cSHong Zhang ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr); 1191c5af039cSHong Zhang } 1192c5af039cSHong Zhang 1193d6ab1dc0SHong Zhang /* zero the current row of Pt*A */ 1194d6ab1dc0SHong Zhang aJ = aj + ai[i]; 1195e745005bSHong Zhang for (k=0; k<anz; k++) aval[aJ[k]] = 0.0; 1196c5af039cSHong Zhang } 1197c5af039cSHong Zhang 1198c5af039cSHong Zhang /* 3) send and recv matrix values coa */ 1199c5af039cSHong Zhang /*------------------------------------*/ 1200c5af039cSHong Zhang buf_ri = merge->buf_ri; 1201c5af039cSHong Zhang buf_rj = merge->buf_rj; 1202c5af039cSHong Zhang len_s = merge->len_s; 1203c5af039cSHong Zhang ierr = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr); 1204c5af039cSHong Zhang ierr = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr); 1205c5af039cSHong Zhang 1206c5af039cSHong Zhang ierr = PetscMalloc2(merge->nsend+1,MPI_Request,&s_waits,size,MPI_Status,&status);CHKERRQ(ierr); 1207c5af039cSHong Zhang for (proc=0,k=0; proc<size; proc++){ 1208c5af039cSHong Zhang if (!len_s[proc]) continue; 1209c5af039cSHong Zhang i = merge->owners_co[proc]; 1210c5af039cSHong Zhang ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr); 1211c5af039cSHong Zhang k++; 1212c5af039cSHong Zhang } 1213c5af039cSHong Zhang if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);} 1214c5af039cSHong Zhang if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);} 1215c5af039cSHong Zhang 1216c5af039cSHong Zhang ierr = PetscFree2(s_waits,status);CHKERRQ(ierr); 1217c5af039cSHong Zhang ierr = PetscFree(r_waits);CHKERRQ(ierr); 1218c5af039cSHong Zhang ierr = PetscFree(coa);CHKERRQ(ierr); 1219c5af039cSHong Zhang 1220c5af039cSHong Zhang /* 4) insert local Cseq and received values into Cmpi */ 1221c5af039cSHong Zhang /*----------------------------------------------------*/ 1222c5af039cSHong Zhang ierr = PetscMalloc3(merge->nrecv,PetscInt**,&buf_ri_k,merge->nrecv,PetscInt*,&nextrow,merge->nrecv,PetscInt*,&nextci);CHKERRQ(ierr); 1223c5af039cSHong Zhang for (k=0; k<merge->nrecv; k++){ 1224c36aecf5SHong Zhang buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 1225c5af039cSHong Zhang nrows = *(buf_ri_k[k]); 1226c5af039cSHong Zhang nextrow[k] = buf_ri_k[k]+1; /* next row number of k-th recved i-structure */ 1227c5af039cSHong Zhang nextci[k] = buf_ri_k[k] + (nrows + 1);/* poins to the next i-structure of k-th recved i-structure */ 1228c5af039cSHong Zhang } 1229c5af039cSHong Zhang 1230c5af039cSHong Zhang for (i=0; i<cm; i++) { 1231c5af039cSHong Zhang row = owners[rank] + i; /* global row index of C_seq */ 1232c5af039cSHong Zhang bj_i = bj + bi[i]; /* col indices of the i-th row of C */ 1233c5af039cSHong Zhang ba_i = ba + bi[i]; 1234c5af039cSHong Zhang bnz = bi[i+1] - bi[i]; 1235c5af039cSHong Zhang /* add received vals into ba */ 1236c5af039cSHong Zhang for (k=0; k<merge->nrecv; k++){ /* k-th received message */ 1237c5af039cSHong Zhang /* i-th row */ 1238c5af039cSHong Zhang if (i == *nextrow[k]) { 1239c5af039cSHong Zhang cnz = *(nextci[k]+1) - *nextci[k]; 1240c5af039cSHong Zhang cj = buf_rj[k] + *(nextci[k]); 1241c5af039cSHong Zhang ca = abuf_r[k] + *(nextci[k]); 1242c5af039cSHong Zhang nextcj = 0; 1243c5af039cSHong Zhang for (j=0; nextcj<cnz; j++){ 1244c5af039cSHong Zhang if (bj_i[j] == cj[nextcj]){ /* bcol == ccol */ 1245c5af039cSHong Zhang ba_i[j] += ca[nextcj++]; 1246c5af039cSHong Zhang } 1247c5af039cSHong Zhang } 1248c5af039cSHong Zhang nextrow[k]++; nextci[k]++; 1249c5af039cSHong Zhang ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr); 1250c5af039cSHong Zhang } 1251c5af039cSHong Zhang } 1252c5af039cSHong Zhang ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr); 1253c5af039cSHong Zhang } 1254c5af039cSHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1255c5af039cSHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1256c5af039cSHong Zhang 1257c5af039cSHong Zhang ierr = PetscFree(ba);CHKERRQ(ierr); 1258c5af039cSHong Zhang ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr); 1259c5af039cSHong Zhang ierr = PetscFree(abuf_r);CHKERRQ(ierr); 1260c5af039cSHong Zhang ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr); 1261e745005bSHong Zhang ierr = PetscFree(aval);CHKERRQ(ierr); 1262f96d37fcSHong Zhang PetscFunctionReturn(0); 1263f96d37fcSHong Zhang } 1264f96d37fcSHong Zhang 1265c5af039cSHong Zhang /* This routine is modified from MatPtAPSymbolic_MPIAIJ_MPIAIJ() */ 1266f96d37fcSHong Zhang #undef __FUNCT__ 1267f96d37fcSHong Zhang #define __FUNCT__ "MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ" 1268f96d37fcSHong Zhang PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(Mat P,Mat A,PetscReal fill,Mat *C) 1269f96d37fcSHong Zhang { 1270f96d37fcSHong Zhang PetscErrorCode ierr; 12714e938277SHong Zhang Mat Cmpi,A_loc,POt,PDt; 1272f96d37fcSHong Zhang Mat_PtAPMPI *ptap; 1273f96d37fcSHong Zhang PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 1274f96d37fcSHong Zhang Mat_MPIAIJ *p=(Mat_MPIAIJ*)P->data,*c; 1275f96d37fcSHong Zhang PetscInt *pdti,*pdtj,*poti,*potj,*ptJ; 1276f96d37fcSHong Zhang PetscInt nnz; 12774e938277SHong Zhang PetscInt *lnk,*owners_co,*coi,*coj,i,k,pnz,row; 1278497f5370SHong Zhang PetscInt am=A->rmap->n,pn=P->cmap->n; 1279f96d37fcSHong Zhang PetscBT lnkbt; 1280f96d37fcSHong Zhang MPI_Comm comm=((PetscObject)A)->comm; 1281f96d37fcSHong Zhang PetscMPIInt size,rank,tagi,tagj,*len_si,*len_s,*len_ri; 1282f96d37fcSHong Zhang PetscInt **buf_rj,**buf_ri,**buf_ri_k; 1283f96d37fcSHong Zhang PetscInt len,proc,*dnz,*onz,*owners; 1284f96d37fcSHong Zhang PetscInt nzi,*bi,*bj; 1285f96d37fcSHong Zhang PetscInt nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci; 1286f96d37fcSHong Zhang MPI_Request *swaits,*rwaits; 1287f96d37fcSHong Zhang MPI_Status *sstatus,rstatus; 1288f96d37fcSHong Zhang Mat_Merge_SeqsToMPI *merge; 1289d6ab1dc0SHong Zhang PetscInt *ai,*aj,*Jptr,anz,*prmap=p->garray,pon,nspacedouble=0,j; 1290f96d37fcSHong Zhang PetscReal afill=1.0,afill_tmp; 1291c36aecf5SHong Zhang PetscInt rstart = P->cmap->rstart,rmax,aN=A->cmap->N,Crmax; 1292f96d37fcSHong Zhang PetscScalar *vals; 12934e938277SHong Zhang Mat_SeqAIJ *a_loc, *pdt,*pot; 1294f96d37fcSHong Zhang 1295f96d37fcSHong Zhang PetscFunctionBegin; 1296c5af039cSHong Zhang /* check if matrix local sizes are compatible */ 1297c5af039cSHong Zhang if (A->rmap->rstart != P->rmap->rstart || A->rmap->rend != P->rmap->rend){ 1298c5af039cSHong Zhang SETERRQ4(comm,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, A (%D, %D) != P (%D,%D)",A->rmap->rstart,A->rmap->rend,P->rmap->rstart,P->rmap->rend); 1299c5af039cSHong Zhang } 1300c5af039cSHong Zhang 1301f96d37fcSHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1302f96d37fcSHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 13031c7d5954SHong Zhang #if defined(DEBUG_MATTrMatMult) 13044e938277SHong Zhang ierr = PetscSynchronizedPrintf(comm,"[%d] TransposeMatMultSymbolic P: %d %d, %d %d; A %d %d, %d %d\n",rank,P->rmap->N,P->cmap->N,P->rmap->n,P->cmap->n,A->rmap->N,aN,A->rmap->n,A->cmap->n); 13051c7d5954SHong Zhang ierr = PetscSynchronizedFlush(comm);CHKERRQ(ierr); 13061c7d5954SHong Zhang #endif 1307f96d37fcSHong Zhang 1308f96d37fcSHong Zhang /* create struct Mat_PtAPMPI and attached it to C later */ 1309f96d37fcSHong Zhang ierr = PetscNew(Mat_PtAPMPI,&ptap);CHKERRQ(ierr); 1310f96d37fcSHong Zhang 1311f96d37fcSHong Zhang /* get A_loc by taking all local rows of A */ 1312f96d37fcSHong Zhang ierr = MatMPIAIJGetLocalMat(A,MAT_INITIAL_MATRIX,&A_loc);CHKERRQ(ierr); 1313c5af039cSHong Zhang ptap->A_loc = A_loc; 13141c7d5954SHong Zhang a_loc = (Mat_SeqAIJ*)(A_loc)->data; 1315d6ab1dc0SHong Zhang ai = a_loc->i; 1316d6ab1dc0SHong Zhang aj = a_loc->j; 1317f96d37fcSHong Zhang 1318f96d37fcSHong Zhang /* determine symbolic Co=(p->B)^T*A - send to others */ 1319f96d37fcSHong Zhang /*----------------------------------------------------*/ 13204e938277SHong Zhang ierr = MatTransposeSymbolic_SeqAIJ(p->A,&PDt);CHKERRQ(ierr); 13214e938277SHong Zhang pdt = (Mat_SeqAIJ*)PDt->data; 13224e938277SHong Zhang pdti = pdt->i; pdtj = pdt->j; 13234e938277SHong Zhang 13244e938277SHong Zhang ierr = MatTransposeSymbolic_SeqAIJ(p->B,&POt);CHKERRQ(ierr); 13254e938277SHong Zhang pot = (Mat_SeqAIJ*)POt->data; 13264e938277SHong Zhang poti = pot->i; potj = pot->j; 1327f96d37fcSHong Zhang 1328f96d37fcSHong Zhang /* then, compute symbolic Co = (p->B)^T*A */ 1329f96d37fcSHong Zhang pon = (p->B)->cmap->n; /* total num of rows to be sent to other processors 1330f96d37fcSHong Zhang >= (num of nonzero rows of C_seq) - pn */ 1331f96d37fcSHong Zhang ierr = PetscMalloc((pon+1)*sizeof(PetscInt),&coi);CHKERRQ(ierr); 1332f96d37fcSHong Zhang coi[0] = 0; 1333f96d37fcSHong Zhang 1334f96d37fcSHong Zhang /* set initial free space to be fill*(nnz(p->B) + nnz(A)) */ 1335d6ab1dc0SHong Zhang nnz = fill*(poti[pon] + ai[am]); 1336f96d37fcSHong Zhang ierr = PetscFreeSpaceGet(nnz,&free_space); 1337f96d37fcSHong Zhang current_space = free_space; 13381c7d5954SHong Zhang #if defined(DEBUG_MATTrMatMult) 1339d6ab1dc0SHong Zhang ierr = PetscSynchronizedPrintf(comm, "[%d] nnz = fill %g *(%d + %d)\n",rank,fill,poti[pon],ai[am]);CHKERRQ(ierr); 1340f96d37fcSHong Zhang ierr = PetscSynchronizedFlush(comm);CHKERRQ(ierr); 13411c7d5954SHong Zhang #endif 1342f96d37fcSHong Zhang /* create and initialize a linked list */ 13434e938277SHong Zhang i = PetscMax(pdt->rmax,pot->rmax); 1344c36aecf5SHong Zhang Crmax = i*a_loc->rmax*size; 13454e938277SHong Zhang if (!Crmax || Crmax > aN) Crmax = aN; 13464e938277SHong Zhang #if defined(DEBUG_MATTrMatMult) 13474e938277SHong Zhang printf("[%d] rmax A_loc %d * max(PD %d, PO %d)=%d, Crmax %d\n",rank,a_loc->rmax,pdt->rmax,pot->rmax,i*a_loc->rmax,Crmax); 13484e938277SHong Zhang #endif 13494e938277SHong Zhang ierr = PetscLLCondensedCreate(Crmax,aN,&lnk,&lnkbt);CHKERRQ(ierr); 1350f96d37fcSHong Zhang 1351f96d37fcSHong Zhang for (i=0; i<pon; i++) { 1352f96d37fcSHong Zhang pnz = poti[i+1] - poti[i]; 1353f96d37fcSHong Zhang ptJ = potj + poti[i]; 1354f96d37fcSHong Zhang for (j=0; j<pnz; j++){ 1355f96d37fcSHong Zhang row = ptJ[j]; /* row of A_loc == col of Pot */ 1356d6ab1dc0SHong Zhang anz = ai[row+1] - ai[row]; 1357d6ab1dc0SHong Zhang Jptr = aj + ai[row]; 1358f96d37fcSHong Zhang /* add non-zero cols of AP into the sorted linked list lnk */ 1359d6ab1dc0SHong Zhang ierr = PetscLLCondensedAddSorted(anz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 1360f96d37fcSHong Zhang } 13614e938277SHong Zhang nnz = lnk[0]; 1362f96d37fcSHong Zhang 1363f96d37fcSHong Zhang /* If free space is not available, double the total space in the list */ 1364f96d37fcSHong Zhang if (current_space->local_remaining<nnz) { 1365f96d37fcSHong Zhang ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 1366f96d37fcSHong Zhang nspacedouble++; 1367f96d37fcSHong Zhang } 1368f96d37fcSHong Zhang 1369f96d37fcSHong Zhang /* Copy data into free space, and zero out denserows */ 13704e938277SHong Zhang ierr = PetscLLCondensedClean(aN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 1371f96d37fcSHong Zhang current_space->array += nnz; 1372f96d37fcSHong Zhang current_space->local_used += nnz; 1373f96d37fcSHong Zhang current_space->local_remaining -= nnz; 1374f96d37fcSHong Zhang coi[i+1] = coi[i] + nnz; 1375f96d37fcSHong Zhang } 1376f96d37fcSHong Zhang 1377f96d37fcSHong Zhang ierr = PetscMalloc((coi[pon]+1)*sizeof(PetscInt),&coj);CHKERRQ(ierr); 1378f96d37fcSHong Zhang ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr); 1379d6ab1dc0SHong Zhang afill_tmp = (PetscReal)coi[pon]/(poti[pon] + ai[am]); 1380f96d37fcSHong Zhang if (afill_tmp > afill) afill = afill_tmp; 1381f96d37fcSHong Zhang 1382f96d37fcSHong Zhang /* send j-array (coj) of Co to other processors */ 1383f96d37fcSHong Zhang /*----------------------------------------------*/ 1384f96d37fcSHong Zhang /* determine row ownership */ 1385f96d37fcSHong Zhang ierr = PetscNew(Mat_Merge_SeqsToMPI,&merge);CHKERRQ(ierr); 1386f96d37fcSHong Zhang ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr); 1387f96d37fcSHong Zhang merge->rowmap->n = pn; 1388f96d37fcSHong Zhang merge->rowmap->bs = 1; 1389f96d37fcSHong Zhang ierr = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr); 1390f96d37fcSHong Zhang owners = merge->rowmap->range; 1391f96d37fcSHong Zhang 1392f96d37fcSHong Zhang /* determine the number of messages to send, their lengths */ 1393f96d37fcSHong Zhang ierr = PetscMalloc(size*sizeof(PetscMPIInt),&len_si);CHKERRQ(ierr); 1394f96d37fcSHong Zhang ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 1395f96d37fcSHong Zhang ierr = PetscMalloc(size*sizeof(PetscMPIInt),&merge->len_s);CHKERRQ(ierr); 1396f96d37fcSHong Zhang len_s = merge->len_s; 1397f96d37fcSHong Zhang merge->nsend = 0; 1398f96d37fcSHong Zhang 1399f96d37fcSHong Zhang ierr = PetscMalloc((size+2)*sizeof(PetscInt),&owners_co);CHKERRQ(ierr); 1400f96d37fcSHong Zhang ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 1401f96d37fcSHong Zhang 1402f96d37fcSHong Zhang proc = 0; 1403f96d37fcSHong Zhang for (i=0; i<pon; i++){ 1404f96d37fcSHong Zhang while (prmap[i] >= owners[proc+1]) proc++; 1405f96d37fcSHong Zhang len_si[proc]++; /* num of rows in Co to be sent to [proc] */ 1406f96d37fcSHong Zhang len_s[proc] += coi[i+1] - coi[i]; 1407f96d37fcSHong Zhang } 1408f96d37fcSHong Zhang 1409f96d37fcSHong Zhang len = 0; /* max length of buf_si[] */ 1410f96d37fcSHong Zhang owners_co[0] = 0; 1411f96d37fcSHong Zhang for (proc=0; proc<size; proc++){ 1412f96d37fcSHong Zhang owners_co[proc+1] = owners_co[proc] + len_si[proc]; 1413f96d37fcSHong Zhang if (len_si[proc]){ 1414f96d37fcSHong Zhang merge->nsend++; 1415f96d37fcSHong Zhang len_si[proc] = 2*(len_si[proc] + 1); 1416f96d37fcSHong Zhang len += len_si[proc]; 1417f96d37fcSHong Zhang } 1418f96d37fcSHong Zhang } 1419f96d37fcSHong Zhang 1420f96d37fcSHong Zhang /* determine the number and length of messages to receive for coi and coj */ 1421f96d37fcSHong Zhang ierr = PetscGatherNumberOfMessages(comm,PETSC_NULL,len_s,&merge->nrecv);CHKERRQ(ierr); 1422f96d37fcSHong Zhang ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr); 1423f96d37fcSHong Zhang 1424f96d37fcSHong Zhang /* post the Irecv and Isend of coj */ 1425f96d37fcSHong Zhang ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr); 1426f96d37fcSHong Zhang ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr); 1427f96d37fcSHong Zhang ierr = PetscMalloc((merge->nsend+1)*sizeof(MPI_Request),&swaits);CHKERRQ(ierr); 1428f96d37fcSHong Zhang for (proc=0, k=0; proc<size; proc++){ 1429f96d37fcSHong Zhang if (!len_s[proc]) continue; 1430f96d37fcSHong Zhang i = owners_co[proc]; 1431f96d37fcSHong Zhang ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr); 1432f96d37fcSHong Zhang k++; 1433f96d37fcSHong Zhang } 1434f96d37fcSHong Zhang 1435f96d37fcSHong Zhang /* receives and sends of coj are complete */ 1436f96d37fcSHong Zhang ierr = PetscMalloc(size*sizeof(MPI_Status),&sstatus);CHKERRQ(ierr); 1437f96d37fcSHong Zhang for (i=0; i<merge->nrecv; i++){ 1438f96d37fcSHong Zhang ierr = MPI_Waitany(merge->nrecv,rwaits,&j,&rstatus);CHKERRQ(ierr); 1439f96d37fcSHong Zhang } 1440f96d37fcSHong Zhang ierr = PetscFree(rwaits);CHKERRQ(ierr); 1441f96d37fcSHong Zhang if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);} 1442f96d37fcSHong Zhang 1443f96d37fcSHong Zhang /* send and recv coi */ 1444f96d37fcSHong Zhang /*-------------------*/ 1445f96d37fcSHong Zhang ierr = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr); 1446f96d37fcSHong Zhang ierr = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr); 1447f96d37fcSHong Zhang ierr = PetscMalloc((len+1)*sizeof(PetscInt),&buf_s);CHKERRQ(ierr); 1448f96d37fcSHong Zhang buf_si = buf_s; /* points to the beginning of k-th msg to be sent */ 1449f96d37fcSHong Zhang for (proc=0,k=0; proc<size; proc++){ 1450f96d37fcSHong Zhang if (!len_s[proc]) continue; 1451f96d37fcSHong Zhang /* form outgoing message for i-structure: 1452f96d37fcSHong Zhang buf_si[0]: nrows to be sent 1453f96d37fcSHong Zhang [1:nrows]: row index (global) 1454f96d37fcSHong Zhang [nrows+1:2*nrows+1]: i-structure index 1455f96d37fcSHong Zhang */ 1456f96d37fcSHong Zhang /*-------------------------------------------*/ 1457f96d37fcSHong Zhang nrows = len_si[proc]/2 - 1; 1458f96d37fcSHong Zhang buf_si_i = buf_si + nrows+1; 1459f96d37fcSHong Zhang buf_si[0] = nrows; 1460f96d37fcSHong Zhang buf_si_i[0] = 0; 1461f96d37fcSHong Zhang nrows = 0; 1462f96d37fcSHong Zhang for (i=owners_co[proc]; i<owners_co[proc+1]; i++){ 1463f96d37fcSHong Zhang nzi = coi[i+1] - coi[i]; 1464f96d37fcSHong Zhang buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */ 1465f96d37fcSHong Zhang buf_si[nrows+1] =prmap[i] -owners[proc]; /* local row index */ 1466f96d37fcSHong Zhang nrows++; 1467f96d37fcSHong Zhang } 1468f96d37fcSHong Zhang ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr); 1469f96d37fcSHong Zhang k++; 1470f96d37fcSHong Zhang buf_si += len_si[proc]; 1471f96d37fcSHong Zhang } 1472f96d37fcSHong Zhang i = merge->nrecv; 1473f96d37fcSHong Zhang while (i--) { 1474f96d37fcSHong Zhang ierr = MPI_Waitany(merge->nrecv,rwaits,&j,&rstatus);CHKERRQ(ierr); 1475f96d37fcSHong Zhang } 1476f96d37fcSHong Zhang ierr = PetscFree(rwaits);CHKERRQ(ierr); 1477f96d37fcSHong Zhang if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);} 1478f96d37fcSHong Zhang ierr = PetscFree(len_si);CHKERRQ(ierr); 1479f96d37fcSHong Zhang ierr = PetscFree(len_ri);CHKERRQ(ierr); 1480f96d37fcSHong Zhang ierr = PetscFree(swaits);CHKERRQ(ierr); 1481f96d37fcSHong Zhang ierr = PetscFree(sstatus);CHKERRQ(ierr); 1482f96d37fcSHong Zhang ierr = PetscFree(buf_s);CHKERRQ(ierr); 1483f96d37fcSHong Zhang 1484f96d37fcSHong Zhang /* compute the local portion of C (mpi mat) */ 1485f96d37fcSHong Zhang /*------------------------------------------*/ 1486f96d37fcSHong Zhang /* allocate bi array and free space for accumulating nonzero column info */ 1487f96d37fcSHong Zhang ierr = PetscMalloc((pn+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr); 1488f96d37fcSHong Zhang bi[0] = 0; 1489f96d37fcSHong Zhang 1490c36aecf5SHong Zhang /* set initial free space to be fill*(nnz(P) + nnz(A)) */ 1491d6ab1dc0SHong Zhang nnz = fill*(pdti[pn] + poti[pon] + ai[am]); 1492f96d37fcSHong Zhang ierr = PetscFreeSpaceGet(nnz,&free_space); 1493f96d37fcSHong Zhang current_space = free_space; 1494f96d37fcSHong Zhang 14951c7d5954SHong Zhang #if defined(DEBUG_MATTrMatMult) 1496c36aecf5SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] nnz=fill %g*(%d + %d + %d)\n",rank,fill,pdti[pn],poti[pon],ai[am]); 14971c7d5954SHong Zhang #endif 1498f96d37fcSHong Zhang 1499f96d37fcSHong Zhang ierr = PetscMalloc3(merge->nrecv,PetscInt**,&buf_ri_k,merge->nrecv,PetscInt*,&nextrow,merge->nrecv,PetscInt*,&nextci);CHKERRQ(ierr); 1500f96d37fcSHong Zhang for (k=0; k<merge->nrecv; k++){ 1501f96d37fcSHong Zhang buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 1502f96d37fcSHong Zhang nrows = *buf_ri_k[k]; 1503f96d37fcSHong Zhang nextrow[k] = buf_ri_k[k] + 1; /* next row number of k-th recved i-structure */ 1504f96d37fcSHong Zhang nextci[k] = buf_ri_k[k] + (nrows + 1);/* poins to the next i-structure of k-th recved i-structure */ 1505f96d37fcSHong Zhang } 1506f96d37fcSHong Zhang 15071c7d5954SHong Zhang ierr = MatPreallocateInitialize(comm,pn,A->cmap->n,dnz,onz);CHKERRQ(ierr); 1508f96d37fcSHong Zhang rmax = 0; 1509f96d37fcSHong Zhang for (i=0; i<pn; i++) { 1510f96d37fcSHong Zhang /* add pdt[i,:]*AP into lnk */ 1511f96d37fcSHong Zhang pnz = pdti[i+1] - pdti[i]; 1512f96d37fcSHong Zhang ptJ = pdtj + pdti[i]; 1513f96d37fcSHong Zhang for (j=0; j<pnz; j++){ 1514f96d37fcSHong Zhang row = ptJ[j]; /* row of AP == col of Pt */ 1515d6ab1dc0SHong Zhang anz = ai[row+1] - ai[row]; 1516d6ab1dc0SHong Zhang Jptr = aj + ai[row]; 1517f96d37fcSHong Zhang /* add non-zero cols of AP into the sorted linked list lnk */ 1518d6ab1dc0SHong Zhang ierr = PetscLLCondensedAddSorted(anz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 1519f96d37fcSHong Zhang } 15204e938277SHong Zhang 1521f96d37fcSHong Zhang /* add received col data into lnk */ 1522f96d37fcSHong Zhang for (k=0; k<merge->nrecv; k++){ /* k-th received message */ 1523f96d37fcSHong Zhang if (i == *nextrow[k]) { /* i-th row */ 1524f96d37fcSHong Zhang nzi = *(nextci[k]+1) - *nextci[k]; 1525f96d37fcSHong Zhang Jptr = buf_rj[k] + *nextci[k]; 15264e938277SHong Zhang ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr); 1527f96d37fcSHong Zhang nextrow[k]++; nextci[k]++; 1528f96d37fcSHong Zhang } 1529f96d37fcSHong Zhang } 15304e938277SHong Zhang nnz = lnk[0]; 1531f96d37fcSHong Zhang 1532f96d37fcSHong Zhang /* if free space is not available, make more free space */ 1533f96d37fcSHong Zhang if (current_space->local_remaining<nnz) { 1534f96d37fcSHong Zhang ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 1535f96d37fcSHong Zhang nspacedouble++; 1536f96d37fcSHong Zhang } 1537f96d37fcSHong Zhang /* copy data into free space, then initialize lnk */ 15384e938277SHong Zhang ierr = PetscLLCondensedClean(aN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 1539f96d37fcSHong Zhang ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr); 1540f96d37fcSHong Zhang current_space->array += nnz; 1541f96d37fcSHong Zhang current_space->local_used += nnz; 1542f96d37fcSHong Zhang current_space->local_remaining -= nnz; 1543f96d37fcSHong Zhang bi[i+1] = bi[i] + nnz; 1544f96d37fcSHong Zhang if (nnz > rmax) rmax = nnz; 1545f96d37fcSHong Zhang } 1546f96d37fcSHong Zhang ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr); 1547f96d37fcSHong Zhang 1548f96d37fcSHong Zhang ierr = PetscMalloc((bi[pn]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr); 1549f96d37fcSHong Zhang ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr); 1550d6ab1dc0SHong Zhang afill_tmp = (PetscReal)bi[pn]/(pdti[pn] + poti[pon] + ai[am]); 1551f96d37fcSHong Zhang if (afill_tmp > afill) afill = afill_tmp; 1552d6ab1dc0SHong Zhang ierr = PetscLLCondensedDestroy(lnk,lnkbt);CHKERRQ(ierr); 15534e938277SHong Zhang ierr = MatDestroy(&POt);CHKERRQ(ierr); 15544e938277SHong Zhang ierr = MatDestroy(&PDt);CHKERRQ(ierr); 1555f96d37fcSHong Zhang 15561c7d5954SHong Zhang /* create symbolic parallel matrix Cmpi - why cannot be assembled in Numeric part */ 15571c7d5954SHong Zhang /*----------------------------------------------------------------------------------*/ 15581c7d5954SHong Zhang ierr = PetscMalloc((rmax+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr); 15591c7d5954SHong Zhang ierr = PetscMemzero(vals,rmax*sizeof(PetscScalar));CHKERRQ(ierr); 1560f96d37fcSHong Zhang 1561f96d37fcSHong Zhang ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr); 1562f96d37fcSHong Zhang ierr = MatSetSizes(Cmpi,pn,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 1563f96d37fcSHong Zhang ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr); 1564f96d37fcSHong Zhang ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr); 1565f96d37fcSHong Zhang ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 1566f96d37fcSHong Zhang ierr = MatSetBlockSize(Cmpi,1);CHKERRQ(ierr); 1567f96d37fcSHong Zhang for (i=0; i<pn; i++){ 1568f96d37fcSHong Zhang row = i + rstart; 1569f96d37fcSHong Zhang nnz = bi[i+1] - bi[i]; 1570f96d37fcSHong Zhang Jptr = bj + bi[i]; 1571f96d37fcSHong Zhang ierr = MatSetValues(Cmpi,1,&row,nnz,Jptr,vals,INSERT_VALUES);CHKERRQ(ierr); 1572f96d37fcSHong Zhang } 1573f96d37fcSHong Zhang ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1574f96d37fcSHong Zhang ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1575f96d37fcSHong Zhang ierr = PetscFree(vals);CHKERRQ(ierr); 1576f96d37fcSHong Zhang 1577f96d37fcSHong Zhang merge->bi = bi; 1578f96d37fcSHong Zhang merge->bj = bj; 1579f96d37fcSHong Zhang merge->coi = coi; 1580f96d37fcSHong Zhang merge->coj = coj; 1581f96d37fcSHong Zhang merge->buf_ri = buf_ri; 1582f96d37fcSHong Zhang merge->buf_rj = buf_rj; 1583f96d37fcSHong Zhang merge->owners_co = owners_co; 1584f96d37fcSHong Zhang merge->destroy = Cmpi->ops->destroy; 1585f96d37fcSHong Zhang merge->duplicate = Cmpi->ops->duplicate; 1586f96d37fcSHong Zhang 1587d6ab1dc0SHong Zhang Cmpi->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ; 1588c5af039cSHong Zhang Cmpi->ops->destroy = MatDestroy_MPIAIJ_PtAP; 1589f96d37fcSHong Zhang 1590f96d37fcSHong Zhang /* attach the supporting struct to Cmpi for reuse */ 1591f96d37fcSHong Zhang c = (Mat_MPIAIJ*)Cmpi->data; 1592f96d37fcSHong Zhang c->ptap = ptap; 1593c5af039cSHong Zhang ptap->api = PETSC_NULL; 1594c5af039cSHong Zhang ptap->apj = PETSC_NULL; 1595f96d37fcSHong Zhang ptap->merge = merge; 1596d6ab1dc0SHong Zhang ptap->rmax = rmax; 1597d6ab1dc0SHong Zhang 1598d6ab1dc0SHong Zhang *C = Cmpi; 1599d6ab1dc0SHong Zhang #if defined(PETSC_USE_INFO) 1600d6ab1dc0SHong Zhang if (bi[pn] != 0) { 1601d6ab1dc0SHong Zhang ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr); 1602d6ab1dc0SHong Zhang ierr = PetscInfo1(Cmpi,"Use MatTransposeMatMult(A,B,MatReuse,%G,&C) for best performance.\n",afill);CHKERRQ(ierr); 1603d6ab1dc0SHong Zhang } else { 1604d6ab1dc0SHong Zhang ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr); 1605d6ab1dc0SHong Zhang } 1606d6ab1dc0SHong Zhang #endif 1607d6ab1dc0SHong Zhang PetscFunctionReturn(0); 1608d6ab1dc0SHong Zhang } 1609d6ab1dc0SHong Zhang 1610d6ab1dc0SHong Zhang #undef __FUNCT__ 1611d6ab1dc0SHong Zhang #define __FUNCT__ "MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_Scalable" 1612d6ab1dc0SHong Zhang PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_Scalable(Mat P,Mat A,Mat C) 1613d6ab1dc0SHong Zhang { 1614d6ab1dc0SHong Zhang PetscErrorCode ierr; 1615d6ab1dc0SHong Zhang Mat_Merge_SeqsToMPI *merge; 1616d6ab1dc0SHong Zhang Mat_MPIAIJ *p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data; 1617d6ab1dc0SHong Zhang Mat_SeqAIJ *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data; 1618d6ab1dc0SHong Zhang Mat_PtAPMPI *ptap; 1619e745005bSHong Zhang PetscInt *adj; 1620e745005bSHong Zhang PetscInt i,j,k,anz,pnz,row,*cj,nexta; 1621e745005bSHong Zhang MatScalar *ada,*ca,valtmp; 1622d6ab1dc0SHong Zhang PetscInt am=A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n; 1623d6ab1dc0SHong Zhang MPI_Comm comm=((PetscObject)C)->comm; 1624d6ab1dc0SHong Zhang PetscMPIInt size,rank,taga,*len_s; 1625d6ab1dc0SHong Zhang PetscInt *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci; 1626d6ab1dc0SHong Zhang PetscInt **buf_ri,**buf_rj; 1627d6ab1dc0SHong Zhang PetscInt cnz=0,*bj_i,*bi,*bj,bnz,nextcj; /* bi,bj,ba: local array of C(mpi mat) */ 1628d6ab1dc0SHong Zhang MPI_Request *s_waits,*r_waits; 1629d6ab1dc0SHong Zhang MPI_Status *status; 1630d6ab1dc0SHong Zhang MatScalar **abuf_r,*ba_i,*pA,*coa,*ba; 1631d6ab1dc0SHong Zhang PetscInt *ai,*aj,*coi,*coj; 1632d6ab1dc0SHong Zhang PetscInt *poJ=po->j,*pdJ=pd->j; 1633d6ab1dc0SHong Zhang Mat A_loc; 1634d6ab1dc0SHong Zhang Mat_SeqAIJ *a_loc; 1635d6ab1dc0SHong Zhang 1636d6ab1dc0SHong Zhang PetscFunctionBegin; 1637d6ab1dc0SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1638d6ab1dc0SHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 1639d6ab1dc0SHong Zhang 1640d6ab1dc0SHong Zhang ptap = c->ptap; 1641d6ab1dc0SHong Zhang merge = ptap->merge; 1642d6ab1dc0SHong Zhang 1643d6ab1dc0SHong Zhang #if defined(DEBUG_MATTrMatMult) 1644d6ab1dc0SHong Zhang ierr = PetscSynchronizedPrintf(comm,"[%d] TransposeMatMultNumeric_Scalable: Crmax %d \n",rank,ptap->rmax); 1645d6ab1dc0SHong Zhang ierr = PetscSynchronizedFlush(comm);CHKERRQ(ierr); 1646d6ab1dc0SHong Zhang #endif 1647d6ab1dc0SHong Zhang 1648e745005bSHong Zhang /* 2) compute numeric C_seq = P_loc^T*A_loc */ 1649e745005bSHong Zhang /*------------------------------------------*/ 1650d6ab1dc0SHong Zhang /* get data from symbolic products */ 1651d6ab1dc0SHong Zhang coi = merge->coi; coj = merge->coj; 1652d6ab1dc0SHong Zhang ierr = PetscMalloc((coi[pon]+1)*sizeof(MatScalar),&coa);CHKERRQ(ierr); 1653d6ab1dc0SHong Zhang ierr = PetscMemzero(coa,coi[pon]*sizeof(MatScalar));CHKERRQ(ierr); 1654d6ab1dc0SHong Zhang 1655d6ab1dc0SHong Zhang bi = merge->bi; bj = merge->bj; 1656d6ab1dc0SHong Zhang owners = merge->rowmap->range; 1657d6ab1dc0SHong Zhang ierr = PetscMalloc((bi[cm]+1)*sizeof(MatScalar),&ba);CHKERRQ(ierr); 1658d6ab1dc0SHong Zhang ierr = PetscMemzero(ba,bi[cm]*sizeof(MatScalar));CHKERRQ(ierr); 1659d6ab1dc0SHong Zhang 1660d6ab1dc0SHong Zhang /* get A_loc by taking all local rows of A */ 1661d6ab1dc0SHong Zhang A_loc = ptap->A_loc; 1662d6ab1dc0SHong Zhang ierr = MatMPIAIJGetLocalMat(A,MAT_REUSE_MATRIX,&A_loc);CHKERRQ(ierr); 1663d6ab1dc0SHong Zhang a_loc = (Mat_SeqAIJ*)(A_loc)->data; 1664d6ab1dc0SHong Zhang ai = a_loc->i; 1665d6ab1dc0SHong Zhang aj = a_loc->j; 1666d6ab1dc0SHong Zhang 1667d6ab1dc0SHong Zhang for (i=0; i<am; i++) { 1668e745005bSHong Zhang /* 2-a) put A[i,:] to dense array aval */ 1669d6ab1dc0SHong Zhang anz = ai[i+1] - ai[i]; 1670d6ab1dc0SHong Zhang adj = aj + ai[i]; 1671d6ab1dc0SHong Zhang ada = a_loc->a + ai[i]; 1672d6ab1dc0SHong Zhang 1673d6ab1dc0SHong Zhang /* 2-b) Compute Cseq = P_loc[i,:]^T*A[i,:] using outer product */ 1674e745005bSHong Zhang /*-------------------------------------------------------------*/ 1675d6ab1dc0SHong Zhang /* put the value into Co=(p->B)^T*A (off-diagonal part, send to others) */ 1676d6ab1dc0SHong Zhang pnz = po->i[i+1] - po->i[i]; 1677d6ab1dc0SHong Zhang poJ = po->j + po->i[i]; 1678d6ab1dc0SHong Zhang pA = po->a + po->i[i]; 1679d6ab1dc0SHong Zhang for (j=0; j<pnz; j++){ 1680d6ab1dc0SHong Zhang row = poJ[j]; 1681d6ab1dc0SHong Zhang cnz = coi[row+1] - coi[row]; 1682d6ab1dc0SHong Zhang cj = coj + coi[row]; 1683d6ab1dc0SHong Zhang ca = coa + coi[row]; 1684e745005bSHong Zhang /* perform sparse axpy */ 1685e745005bSHong Zhang nexta = 0; 1686d6ab1dc0SHong Zhang valtmp = pA[j]; 1687e745005bSHong Zhang for (k=0; nexta<anz; k++) { 1688e745005bSHong Zhang if (cj[k] == adj[nexta]){ 1689e745005bSHong Zhang ca[k] += valtmp*ada[nexta]; 1690e745005bSHong Zhang nexta++; 1691d6ab1dc0SHong Zhang } 1692e745005bSHong Zhang } 1693e745005bSHong Zhang ierr = PetscLogFlops(2.0*anz);CHKERRQ(ierr); 1694d6ab1dc0SHong Zhang } 1695d6ab1dc0SHong Zhang 1696d6ab1dc0SHong Zhang /* put the value into Cd (diagonal part) */ 1697d6ab1dc0SHong Zhang pnz = pd->i[i+1] - pd->i[i]; 1698d6ab1dc0SHong Zhang pdJ = pd->j + pd->i[i]; 1699d6ab1dc0SHong Zhang pA = pd->a + pd->i[i]; 1700d6ab1dc0SHong Zhang for (j=0; j<pnz; j++){ 1701d6ab1dc0SHong Zhang row = pdJ[j]; 1702d6ab1dc0SHong Zhang cnz = bi[row+1] - bi[row]; 1703d6ab1dc0SHong Zhang cj = bj + bi[row]; 1704d6ab1dc0SHong Zhang ca = ba + bi[row]; 1705e745005bSHong Zhang /* perform sparse axpy */ 1706e745005bSHong Zhang nexta = 0; 1707d6ab1dc0SHong Zhang valtmp = pA[j]; 1708e745005bSHong Zhang for (k=0; nexta<anz; k++) { 1709e745005bSHong Zhang if (cj[k] == adj[nexta]){ 1710e745005bSHong Zhang ca[k] += valtmp*ada[nexta]; 1711e745005bSHong Zhang nexta++; 1712d6ab1dc0SHong Zhang } 1713e745005bSHong Zhang } 1714e745005bSHong Zhang ierr = PetscLogFlops(2.0*anz);CHKERRQ(ierr); 1715d6ab1dc0SHong Zhang } 1716d6ab1dc0SHong Zhang 1717d6ab1dc0SHong Zhang } 1718d6ab1dc0SHong Zhang 1719d6ab1dc0SHong Zhang /* 3) send and recv matrix values coa */ 1720d6ab1dc0SHong Zhang /*------------------------------------*/ 1721d6ab1dc0SHong Zhang buf_ri = merge->buf_ri; 1722d6ab1dc0SHong Zhang buf_rj = merge->buf_rj; 1723d6ab1dc0SHong Zhang len_s = merge->len_s; 1724d6ab1dc0SHong Zhang ierr = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr); 1725d6ab1dc0SHong Zhang ierr = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr); 1726d6ab1dc0SHong Zhang 1727d6ab1dc0SHong Zhang ierr = PetscMalloc2(merge->nsend+1,MPI_Request,&s_waits,size,MPI_Status,&status);CHKERRQ(ierr); 1728d6ab1dc0SHong Zhang for (proc=0,k=0; proc<size; proc++){ 1729d6ab1dc0SHong Zhang if (!len_s[proc]) continue; 1730d6ab1dc0SHong Zhang i = merge->owners_co[proc]; 1731d6ab1dc0SHong Zhang ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr); 1732d6ab1dc0SHong Zhang k++; 1733d6ab1dc0SHong Zhang } 1734d6ab1dc0SHong Zhang if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);} 1735d6ab1dc0SHong Zhang if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);} 1736d6ab1dc0SHong Zhang 1737d6ab1dc0SHong Zhang ierr = PetscFree2(s_waits,status);CHKERRQ(ierr); 1738d6ab1dc0SHong Zhang ierr = PetscFree(r_waits);CHKERRQ(ierr); 1739d6ab1dc0SHong Zhang ierr = PetscFree(coa);CHKERRQ(ierr); 1740d6ab1dc0SHong Zhang 1741d6ab1dc0SHong Zhang /* 4) insert local Cseq and received values into Cmpi */ 1742d6ab1dc0SHong Zhang /*----------------------------------------------------*/ 1743d6ab1dc0SHong Zhang ierr = PetscMalloc3(merge->nrecv,PetscInt**,&buf_ri_k,merge->nrecv,PetscInt*,&nextrow,merge->nrecv,PetscInt*,&nextci);CHKERRQ(ierr); 1744d6ab1dc0SHong Zhang for (k=0; k<merge->nrecv; k++){ 1745e745005bSHong Zhang buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 1746d6ab1dc0SHong Zhang nrows = *(buf_ri_k[k]); 1747d6ab1dc0SHong Zhang nextrow[k] = buf_ri_k[k]+1; /* next row number of k-th recved i-structure */ 1748d6ab1dc0SHong Zhang nextci[k] = buf_ri_k[k] + (nrows + 1);/* poins to the next i-structure of k-th recved i-structure */ 1749d6ab1dc0SHong Zhang } 1750d6ab1dc0SHong Zhang 1751d6ab1dc0SHong Zhang for (i=0; i<cm; i++) { 1752d6ab1dc0SHong Zhang row = owners[rank] + i; /* global row index of C_seq */ 1753d6ab1dc0SHong Zhang bj_i = bj + bi[i]; /* col indices of the i-th row of C */ 1754d6ab1dc0SHong Zhang ba_i = ba + bi[i]; 1755d6ab1dc0SHong Zhang bnz = bi[i+1] - bi[i]; 1756d6ab1dc0SHong Zhang /* add received vals into ba */ 1757d6ab1dc0SHong Zhang for (k=0; k<merge->nrecv; k++){ /* k-th received message */ 1758d6ab1dc0SHong Zhang /* i-th row */ 1759d6ab1dc0SHong Zhang if (i == *nextrow[k]) { 1760d6ab1dc0SHong Zhang cnz = *(nextci[k]+1) - *nextci[k]; 1761d6ab1dc0SHong Zhang cj = buf_rj[k] + *(nextci[k]); 1762d6ab1dc0SHong Zhang ca = abuf_r[k] + *(nextci[k]); 1763d6ab1dc0SHong Zhang nextcj = 0; 1764d6ab1dc0SHong Zhang for (j=0; nextcj<cnz; j++){ 1765d6ab1dc0SHong Zhang if (bj_i[j] == cj[nextcj]){ /* bcol == ccol */ 1766d6ab1dc0SHong Zhang ba_i[j] += ca[nextcj++]; 1767d6ab1dc0SHong Zhang } 1768d6ab1dc0SHong Zhang } 1769d6ab1dc0SHong Zhang nextrow[k]++; nextci[k]++; 1770d6ab1dc0SHong Zhang ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr); 1771d6ab1dc0SHong Zhang } 1772d6ab1dc0SHong Zhang } 1773d6ab1dc0SHong Zhang ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr); 1774d6ab1dc0SHong Zhang } 1775d6ab1dc0SHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1776d6ab1dc0SHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1777d6ab1dc0SHong Zhang 1778d6ab1dc0SHong Zhang ierr = PetscFree(ba);CHKERRQ(ierr); 1779d6ab1dc0SHong Zhang ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr); 1780d6ab1dc0SHong Zhang ierr = PetscFree(abuf_r);CHKERRQ(ierr); 1781d6ab1dc0SHong Zhang ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr); 1782d6ab1dc0SHong Zhang PetscFunctionReturn(0); 1783d6ab1dc0SHong Zhang } 1784d6ab1dc0SHong Zhang 1785d6ab1dc0SHong Zhang /* This routine is modified from MatPtAPSymbolic_MPIAIJ_MPIAIJ() */ 1786d6ab1dc0SHong Zhang #undef __FUNCT__ 1787d6ab1dc0SHong Zhang #define __FUNCT__ "MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable" 1788d6ab1dc0SHong Zhang PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable(Mat P,Mat A,PetscReal fill,Mat *C) 1789d6ab1dc0SHong Zhang { 1790d6ab1dc0SHong Zhang PetscErrorCode ierr; 1791d6ab1dc0SHong Zhang Mat Cmpi,A_loc,POt,PDt; 1792d6ab1dc0SHong Zhang Mat_PtAPMPI *ptap; 1793d6ab1dc0SHong Zhang PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 1794d6ab1dc0SHong Zhang Mat_MPIAIJ *p=(Mat_MPIAIJ*)P->data,*c; 1795d6ab1dc0SHong Zhang PetscInt *pdti,*pdtj,*poti,*potj,*ptJ; 1796d6ab1dc0SHong Zhang PetscInt nnz; 1797d6ab1dc0SHong Zhang PetscInt *lnk,*owners_co,*coi,*coj,i,k,pnz,row; 1798d6ab1dc0SHong Zhang PetscInt am=A->rmap->n,pn=P->cmap->n; 1799d6ab1dc0SHong Zhang MPI_Comm comm=((PetscObject)A)->comm; 1800d6ab1dc0SHong Zhang PetscMPIInt size,rank,tagi,tagj,*len_si,*len_s,*len_ri; 1801d6ab1dc0SHong Zhang PetscInt **buf_rj,**buf_ri,**buf_ri_k; 1802d6ab1dc0SHong Zhang PetscInt len,proc,*dnz,*onz,*owners; 1803d6ab1dc0SHong Zhang PetscInt nzi,*bi,*bj; 1804d6ab1dc0SHong Zhang PetscInt nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci; 1805d6ab1dc0SHong Zhang MPI_Request *swaits,*rwaits; 1806d6ab1dc0SHong Zhang MPI_Status *sstatus,rstatus; 1807d6ab1dc0SHong Zhang Mat_Merge_SeqsToMPI *merge; 1808d6ab1dc0SHong Zhang PetscInt *ai,*aj,*Jptr,anz,*prmap=p->garray,pon,nspacedouble=0,j; 1809d6ab1dc0SHong Zhang PetscReal afill=1.0,afill_tmp; 1810c36aecf5SHong Zhang PetscInt rstart = P->cmap->rstart,rmax,aN=A->cmap->N,Crmax; 1811d6ab1dc0SHong Zhang PetscScalar *vals; 1812d6ab1dc0SHong Zhang Mat_SeqAIJ *a_loc, *pdt,*pot; 1813d6ab1dc0SHong Zhang 1814d6ab1dc0SHong Zhang PetscFunctionBegin; 1815d6ab1dc0SHong Zhang /* check if matrix local sizes are compatible */ 1816d6ab1dc0SHong Zhang if (A->rmap->rstart != P->rmap->rstart || A->rmap->rend != P->rmap->rend){ 1817d6ab1dc0SHong Zhang SETERRQ4(comm,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, A (%D, %D) != P (%D,%D)",A->rmap->rstart,A->rmap->rend,P->rmap->rstart,P->rmap->rend); 1818d6ab1dc0SHong Zhang } 1819d6ab1dc0SHong Zhang 1820d6ab1dc0SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1821d6ab1dc0SHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 1822d6ab1dc0SHong Zhang #if defined(DEBUG_MATTrMatMult) 1823d6ab1dc0SHong Zhang ierr = PetscSynchronizedPrintf(comm,"[%d] TransposeMatMultSymbolic_Scalable P: %d %d, %d %d; A %d %d, %d %d\n",rank,P->rmap->N,P->cmap->N,P->rmap->n,P->cmap->n,A->rmap->N,aN,A->rmap->n,A->cmap->n); 1824d6ab1dc0SHong Zhang ierr = PetscSynchronizedFlush(comm);CHKERRQ(ierr); 1825d6ab1dc0SHong Zhang #endif 1826d6ab1dc0SHong Zhang 1827d6ab1dc0SHong Zhang /* create struct Mat_PtAPMPI and attached it to C later */ 1828d6ab1dc0SHong Zhang ierr = PetscNew(Mat_PtAPMPI,&ptap);CHKERRQ(ierr); 1829d6ab1dc0SHong Zhang 1830d6ab1dc0SHong Zhang /* get A_loc by taking all local rows of A */ 1831d6ab1dc0SHong Zhang ierr = MatMPIAIJGetLocalMat(A,MAT_INITIAL_MATRIX,&A_loc);CHKERRQ(ierr); 1832d6ab1dc0SHong Zhang ptap->A_loc = A_loc; 1833d6ab1dc0SHong Zhang a_loc = (Mat_SeqAIJ*)(A_loc)->data; 1834d6ab1dc0SHong Zhang ai = a_loc->i; 1835d6ab1dc0SHong Zhang aj = a_loc->j; 1836d6ab1dc0SHong Zhang 1837d6ab1dc0SHong Zhang /* determine symbolic Co=(p->B)^T*A - send to others */ 1838d6ab1dc0SHong Zhang /*----------------------------------------------------*/ 1839d6ab1dc0SHong Zhang ierr = MatTransposeSymbolic_SeqAIJ(p->A,&PDt);CHKERRQ(ierr); 1840d6ab1dc0SHong Zhang pdt = (Mat_SeqAIJ*)PDt->data; 1841d6ab1dc0SHong Zhang pdti = pdt->i; pdtj = pdt->j; 1842d6ab1dc0SHong Zhang 1843d6ab1dc0SHong Zhang ierr = MatTransposeSymbolic_SeqAIJ(p->B,&POt);CHKERRQ(ierr); 1844d6ab1dc0SHong Zhang pot = (Mat_SeqAIJ*)POt->data; 1845d6ab1dc0SHong Zhang poti = pot->i; potj = pot->j; 1846d6ab1dc0SHong Zhang 1847d6ab1dc0SHong Zhang /* then, compute symbolic Co = (p->B)^T*A */ 1848d6ab1dc0SHong Zhang pon = (p->B)->cmap->n; /* total num of rows to be sent to other processors 1849d6ab1dc0SHong Zhang >= (num of nonzero rows of C_seq) - pn */ 1850d6ab1dc0SHong Zhang ierr = PetscMalloc((pon+1)*sizeof(PetscInt),&coi);CHKERRQ(ierr); 1851d6ab1dc0SHong Zhang coi[0] = 0; 1852d6ab1dc0SHong Zhang 1853d6ab1dc0SHong Zhang /* set initial free space to be fill*(nnz(p->B) + nnz(A)) */ 1854d6ab1dc0SHong Zhang nnz = fill*(poti[pon] + ai[am]); 1855d6ab1dc0SHong Zhang ierr = PetscFreeSpaceGet(nnz,&free_space); 1856d6ab1dc0SHong Zhang current_space = free_space; 1857d6ab1dc0SHong Zhang #if defined(DEBUG_MATTrMatMult) 1858d6ab1dc0SHong Zhang ierr = PetscSynchronizedPrintf(comm, "[%d] nnz = fill %g *(%d + %d)\n",rank,fill,poti[pon],ai[am]);CHKERRQ(ierr); 1859d6ab1dc0SHong Zhang ierr = PetscSynchronizedFlush(comm);CHKERRQ(ierr); 1860d6ab1dc0SHong Zhang #endif 1861d6ab1dc0SHong Zhang /* create and initialize a linked list */ 1862d6ab1dc0SHong Zhang i = PetscMax(pdt->rmax,pot->rmax); 1863c36aecf5SHong Zhang Crmax = i*a_loc->rmax*size; /* non-scalable! */ 1864d6ab1dc0SHong Zhang if (!Crmax || Crmax > aN) Crmax = aN; 1865d6ab1dc0SHong Zhang #if defined(DEBUG_MATTrMatMult) 1866d6ab1dc0SHong Zhang printf("[%d] rmax A_loc %d * max(PD %d, PO %d)=%d, Crmax %d\n",rank,a_loc->rmax,pdt->rmax,pot->rmax,i*a_loc->rmax,Crmax); 1867d6ab1dc0SHong Zhang #endif 1868d6ab1dc0SHong Zhang ierr = PetscLLCondensedCreate_Scalable(Crmax,&lnk);CHKERRQ(ierr); 1869d6ab1dc0SHong Zhang 1870d6ab1dc0SHong Zhang for (i=0; i<pon; i++) { 1871d6ab1dc0SHong Zhang nnz = 0; 1872d6ab1dc0SHong Zhang pnz = poti[i+1] - poti[i]; 1873d6ab1dc0SHong Zhang ptJ = potj + poti[i]; 1874d6ab1dc0SHong Zhang for (j=0; j<pnz; j++){ 1875d6ab1dc0SHong Zhang row = ptJ[j]; /* row of A_loc == col of Pot */ 1876d6ab1dc0SHong Zhang anz = ai[row+1] - ai[row]; 1877d6ab1dc0SHong Zhang Jptr = aj + ai[row]; 1878d6ab1dc0SHong Zhang /* add non-zero cols of AP into the sorted linked list lnk */ 1879d6ab1dc0SHong Zhang ierr = PetscLLCondensedAddSorted_Scalable(anz,Jptr,lnk);CHKERRQ(ierr); 1880d6ab1dc0SHong Zhang } 1881d6ab1dc0SHong Zhang nnz = lnk[0]; 1882d6ab1dc0SHong Zhang 1883d6ab1dc0SHong Zhang /* If free space is not available, double the total space in the list */ 1884d6ab1dc0SHong Zhang if (current_space->local_remaining<nnz) { 1885d6ab1dc0SHong Zhang ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 1886d6ab1dc0SHong Zhang nspacedouble++; 1887d6ab1dc0SHong Zhang } 1888d6ab1dc0SHong Zhang 1889d6ab1dc0SHong Zhang /* Copy data into free space, and zero out denserows */ 1890d6ab1dc0SHong Zhang ierr = PetscLLCondensedClean_Scalable(nnz,current_space->array,lnk);CHKERRQ(ierr); 1891d6ab1dc0SHong Zhang current_space->array += nnz; 1892d6ab1dc0SHong Zhang current_space->local_used += nnz; 1893d6ab1dc0SHong Zhang current_space->local_remaining -= nnz; 1894d6ab1dc0SHong Zhang coi[i+1] = coi[i] + nnz; 1895d6ab1dc0SHong Zhang } 1896d6ab1dc0SHong Zhang 1897d6ab1dc0SHong Zhang ierr = PetscMalloc((coi[pon]+1)*sizeof(PetscInt),&coj);CHKERRQ(ierr); 1898d6ab1dc0SHong Zhang ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr); 1899d6ab1dc0SHong Zhang afill_tmp = (PetscReal)coi[pon]/(poti[pon] + ai[am]); 1900d6ab1dc0SHong Zhang if (afill_tmp > afill) afill = afill_tmp; 1901d6ab1dc0SHong Zhang 1902d6ab1dc0SHong Zhang /* send j-array (coj) of Co to other processors */ 1903d6ab1dc0SHong Zhang /*----------------------------------------------*/ 1904d6ab1dc0SHong Zhang /* determine row ownership */ 1905d6ab1dc0SHong Zhang ierr = PetscNew(Mat_Merge_SeqsToMPI,&merge);CHKERRQ(ierr); 1906d6ab1dc0SHong Zhang ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr); 1907d6ab1dc0SHong Zhang merge->rowmap->n = pn; 1908d6ab1dc0SHong Zhang merge->rowmap->bs = 1; 1909d6ab1dc0SHong Zhang ierr = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr); 1910d6ab1dc0SHong Zhang owners = merge->rowmap->range; 1911d6ab1dc0SHong Zhang 1912d6ab1dc0SHong Zhang /* determine the number of messages to send, their lengths */ 1913d6ab1dc0SHong Zhang ierr = PetscMalloc(size*sizeof(PetscMPIInt),&len_si);CHKERRQ(ierr); 1914d6ab1dc0SHong Zhang ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 1915d6ab1dc0SHong Zhang ierr = PetscMalloc(size*sizeof(PetscMPIInt),&merge->len_s);CHKERRQ(ierr); 1916d6ab1dc0SHong Zhang len_s = merge->len_s; 1917d6ab1dc0SHong Zhang merge->nsend = 0; 1918d6ab1dc0SHong Zhang 1919d6ab1dc0SHong Zhang ierr = PetscMalloc((size+2)*sizeof(PetscInt),&owners_co);CHKERRQ(ierr); 1920d6ab1dc0SHong Zhang ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 1921d6ab1dc0SHong Zhang 1922d6ab1dc0SHong Zhang proc = 0; 1923d6ab1dc0SHong Zhang for (i=0; i<pon; i++){ 1924d6ab1dc0SHong Zhang while (prmap[i] >= owners[proc+1]) proc++; 1925d6ab1dc0SHong Zhang len_si[proc]++; /* num of rows in Co to be sent to [proc] */ 1926d6ab1dc0SHong Zhang len_s[proc] += coi[i+1] - coi[i]; 1927d6ab1dc0SHong Zhang } 1928d6ab1dc0SHong Zhang 1929d6ab1dc0SHong Zhang len = 0; /* max length of buf_si[] */ 1930d6ab1dc0SHong Zhang owners_co[0] = 0; 1931d6ab1dc0SHong Zhang for (proc=0; proc<size; proc++){ 1932d6ab1dc0SHong Zhang owners_co[proc+1] = owners_co[proc] + len_si[proc]; 1933d6ab1dc0SHong Zhang if (len_si[proc]){ 1934d6ab1dc0SHong Zhang merge->nsend++; 1935d6ab1dc0SHong Zhang len_si[proc] = 2*(len_si[proc] + 1); 1936d6ab1dc0SHong Zhang len += len_si[proc]; 1937d6ab1dc0SHong Zhang } 1938d6ab1dc0SHong Zhang } 1939d6ab1dc0SHong Zhang 1940d6ab1dc0SHong Zhang /* determine the number and length of messages to receive for coi and coj */ 1941d6ab1dc0SHong Zhang ierr = PetscGatherNumberOfMessages(comm,PETSC_NULL,len_s,&merge->nrecv);CHKERRQ(ierr); 1942d6ab1dc0SHong Zhang ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr); 1943d6ab1dc0SHong Zhang 1944d6ab1dc0SHong Zhang /* post the Irecv and Isend of coj */ 1945d6ab1dc0SHong Zhang ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr); 1946d6ab1dc0SHong Zhang ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr); 1947d6ab1dc0SHong Zhang ierr = PetscMalloc((merge->nsend+1)*sizeof(MPI_Request),&swaits);CHKERRQ(ierr); 1948d6ab1dc0SHong Zhang for (proc=0, k=0; proc<size; proc++){ 1949d6ab1dc0SHong Zhang if (!len_s[proc]) continue; 1950d6ab1dc0SHong Zhang i = owners_co[proc]; 1951d6ab1dc0SHong Zhang ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr); 1952d6ab1dc0SHong Zhang k++; 1953d6ab1dc0SHong Zhang } 1954d6ab1dc0SHong Zhang 1955d6ab1dc0SHong Zhang /* receives and sends of coj are complete */ 1956d6ab1dc0SHong Zhang ierr = PetscMalloc(size*sizeof(MPI_Status),&sstatus);CHKERRQ(ierr); 1957d6ab1dc0SHong Zhang for (i=0; i<merge->nrecv; i++){ 1958d6ab1dc0SHong Zhang ierr = MPI_Waitany(merge->nrecv,rwaits,&j,&rstatus);CHKERRQ(ierr); 1959d6ab1dc0SHong Zhang } 1960d6ab1dc0SHong Zhang ierr = PetscFree(rwaits);CHKERRQ(ierr); 1961d6ab1dc0SHong Zhang if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);} 1962d6ab1dc0SHong Zhang 1963d6ab1dc0SHong Zhang /* send and recv coi */ 1964d6ab1dc0SHong Zhang /*-------------------*/ 1965d6ab1dc0SHong Zhang ierr = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr); 1966d6ab1dc0SHong Zhang ierr = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr); 1967d6ab1dc0SHong Zhang ierr = PetscMalloc((len+1)*sizeof(PetscInt),&buf_s);CHKERRQ(ierr); 1968d6ab1dc0SHong Zhang buf_si = buf_s; /* points to the beginning of k-th msg to be sent */ 1969d6ab1dc0SHong Zhang for (proc=0,k=0; proc<size; proc++){ 1970d6ab1dc0SHong Zhang if (!len_s[proc]) continue; 1971d6ab1dc0SHong Zhang /* form outgoing message for i-structure: 1972d6ab1dc0SHong Zhang buf_si[0]: nrows to be sent 1973d6ab1dc0SHong Zhang [1:nrows]: row index (global) 1974d6ab1dc0SHong Zhang [nrows+1:2*nrows+1]: i-structure index 1975d6ab1dc0SHong Zhang */ 1976d6ab1dc0SHong Zhang /*-------------------------------------------*/ 1977d6ab1dc0SHong Zhang nrows = len_si[proc]/2 - 1; 1978d6ab1dc0SHong Zhang buf_si_i = buf_si + nrows+1; 1979d6ab1dc0SHong Zhang buf_si[0] = nrows; 1980d6ab1dc0SHong Zhang buf_si_i[0] = 0; 1981d6ab1dc0SHong Zhang nrows = 0; 1982d6ab1dc0SHong Zhang for (i=owners_co[proc]; i<owners_co[proc+1]; i++){ 1983d6ab1dc0SHong Zhang nzi = coi[i+1] - coi[i]; 1984d6ab1dc0SHong Zhang buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */ 1985d6ab1dc0SHong Zhang buf_si[nrows+1] =prmap[i] -owners[proc]; /* local row index */ 1986d6ab1dc0SHong Zhang nrows++; 1987d6ab1dc0SHong Zhang } 1988d6ab1dc0SHong Zhang ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr); 1989d6ab1dc0SHong Zhang k++; 1990d6ab1dc0SHong Zhang buf_si += len_si[proc]; 1991d6ab1dc0SHong Zhang } 1992d6ab1dc0SHong Zhang i = merge->nrecv; 1993d6ab1dc0SHong Zhang while (i--) { 1994d6ab1dc0SHong Zhang ierr = MPI_Waitany(merge->nrecv,rwaits,&j,&rstatus);CHKERRQ(ierr); 1995d6ab1dc0SHong Zhang } 1996d6ab1dc0SHong Zhang ierr = PetscFree(rwaits);CHKERRQ(ierr); 1997d6ab1dc0SHong Zhang if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);} 1998d6ab1dc0SHong Zhang ierr = PetscFree(len_si);CHKERRQ(ierr); 1999d6ab1dc0SHong Zhang ierr = PetscFree(len_ri);CHKERRQ(ierr); 2000d6ab1dc0SHong Zhang ierr = PetscFree(swaits);CHKERRQ(ierr); 2001d6ab1dc0SHong Zhang ierr = PetscFree(sstatus);CHKERRQ(ierr); 2002d6ab1dc0SHong Zhang ierr = PetscFree(buf_s);CHKERRQ(ierr); 2003d6ab1dc0SHong Zhang 2004d6ab1dc0SHong Zhang /* compute the local portion of C (mpi mat) */ 2005d6ab1dc0SHong Zhang /*------------------------------------------*/ 2006d6ab1dc0SHong Zhang /* allocate bi array and free space for accumulating nonzero column info */ 2007d6ab1dc0SHong Zhang ierr = PetscMalloc((pn+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr); 2008d6ab1dc0SHong Zhang bi[0] = 0; 2009d6ab1dc0SHong Zhang 2010d6ab1dc0SHong Zhang /* set initial free space to be fill*(nnz(P) + nnz(AP)) */ 2011d6ab1dc0SHong Zhang nnz = fill*(pdti[pn] + poti[pon] + ai[am]); 2012d6ab1dc0SHong Zhang ierr = PetscFreeSpaceGet(nnz,&free_space); 2013d6ab1dc0SHong Zhang current_space = free_space; 2014d6ab1dc0SHong Zhang 2015d6ab1dc0SHong Zhang #if defined(DEBUG_MATTrMatMult) 2016d6ab1dc0SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] nnz=%d + %d + %d\n",rank,pdti[pn],poti[pon],ai[am]); 2017d6ab1dc0SHong Zhang #endif 2018d6ab1dc0SHong Zhang 2019d6ab1dc0SHong Zhang ierr = PetscMalloc3(merge->nrecv,PetscInt**,&buf_ri_k,merge->nrecv,PetscInt*,&nextrow,merge->nrecv,PetscInt*,&nextci);CHKERRQ(ierr); 2020d6ab1dc0SHong Zhang for (k=0; k<merge->nrecv; k++){ 2021d6ab1dc0SHong Zhang buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 2022d6ab1dc0SHong Zhang nrows = *buf_ri_k[k]; 2023d6ab1dc0SHong Zhang nextrow[k] = buf_ri_k[k] + 1; /* next row number of k-th recved i-structure */ 2024d6ab1dc0SHong Zhang nextci[k] = buf_ri_k[k] + (nrows + 1);/* poins to the next i-structure of k-th recved i-structure */ 2025d6ab1dc0SHong Zhang } 2026d6ab1dc0SHong Zhang 2027d6ab1dc0SHong Zhang ierr = MatPreallocateInitialize(comm,pn,A->cmap->n,dnz,onz);CHKERRQ(ierr); 2028d6ab1dc0SHong Zhang rmax = 0; 2029d6ab1dc0SHong Zhang for (i=0; i<pn; i++) { 2030d6ab1dc0SHong Zhang /* add pdt[i,:]*AP into lnk */ 2031d6ab1dc0SHong Zhang pnz = pdti[i+1] - pdti[i]; 2032d6ab1dc0SHong Zhang ptJ = pdtj + pdti[i]; 2033d6ab1dc0SHong Zhang for (j=0; j<pnz; j++){ 2034d6ab1dc0SHong Zhang row = ptJ[j]; /* row of AP == col of Pt */ 2035d6ab1dc0SHong Zhang anz = ai[row+1] - ai[row]; 2036d6ab1dc0SHong Zhang Jptr = aj + ai[row]; 2037d6ab1dc0SHong Zhang /* add non-zero cols of AP into the sorted linked list lnk */ 2038d6ab1dc0SHong Zhang ierr = PetscLLCondensedAddSorted_Scalable(anz,Jptr,lnk);CHKERRQ(ierr); 2039d6ab1dc0SHong Zhang } 2040d6ab1dc0SHong Zhang 2041d6ab1dc0SHong Zhang /* add received col data into lnk */ 2042d6ab1dc0SHong Zhang for (k=0; k<merge->nrecv; k++){ /* k-th received message */ 2043d6ab1dc0SHong Zhang if (i == *nextrow[k]) { /* i-th row */ 2044d6ab1dc0SHong Zhang nzi = *(nextci[k]+1) - *nextci[k]; 2045d6ab1dc0SHong Zhang Jptr = buf_rj[k] + *nextci[k]; 2046d6ab1dc0SHong Zhang ierr = PetscLLCondensedAddSorted_Scalable(nzi,Jptr,lnk);CHKERRQ(ierr); 2047d6ab1dc0SHong Zhang nextrow[k]++; nextci[k]++; 2048d6ab1dc0SHong Zhang } 2049d6ab1dc0SHong Zhang } 2050d6ab1dc0SHong Zhang nnz = lnk[0]; 2051d6ab1dc0SHong Zhang 2052d6ab1dc0SHong Zhang /* if free space is not available, make more free space */ 2053d6ab1dc0SHong Zhang if (current_space->local_remaining<nnz) { 2054d6ab1dc0SHong Zhang ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 2055d6ab1dc0SHong Zhang nspacedouble++; 2056d6ab1dc0SHong Zhang } 2057d6ab1dc0SHong Zhang /* copy data into free space, then initialize lnk */ 2058d6ab1dc0SHong Zhang ierr = PetscLLCondensedClean_Scalable(nnz,current_space->array,lnk);CHKERRQ(ierr); 2059d6ab1dc0SHong Zhang ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr); 2060d6ab1dc0SHong Zhang current_space->array += nnz; 2061d6ab1dc0SHong Zhang current_space->local_used += nnz; 2062d6ab1dc0SHong Zhang current_space->local_remaining -= nnz; 2063d6ab1dc0SHong Zhang bi[i+1] = bi[i] + nnz; 2064d6ab1dc0SHong Zhang if (nnz > rmax) rmax = nnz; 2065d6ab1dc0SHong Zhang } 2066d6ab1dc0SHong Zhang ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr); 2067d6ab1dc0SHong Zhang 2068d6ab1dc0SHong Zhang ierr = PetscMalloc((bi[pn]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr); 2069d6ab1dc0SHong Zhang ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr); 2070d6ab1dc0SHong Zhang afill_tmp = (PetscReal)bi[pn]/(pdti[pn] + poti[pon] + ai[am]); 2071d6ab1dc0SHong Zhang if (afill_tmp > afill) afill = afill_tmp; 2072d6ab1dc0SHong Zhang ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr); 2073d6ab1dc0SHong Zhang ierr = MatDestroy(&POt);CHKERRQ(ierr); 2074d6ab1dc0SHong Zhang ierr = MatDestroy(&PDt);CHKERRQ(ierr); 2075d6ab1dc0SHong Zhang 2076d6ab1dc0SHong Zhang /* create symbolic parallel matrix Cmpi - why cannot be assembled in Numeric part */ 2077d6ab1dc0SHong Zhang /*----------------------------------------------------------------------------------*/ 2078d6ab1dc0SHong Zhang ierr = PetscMalloc((rmax+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr); 2079d6ab1dc0SHong Zhang ierr = PetscMemzero(vals,rmax*sizeof(PetscScalar));CHKERRQ(ierr); 2080d6ab1dc0SHong Zhang 2081d6ab1dc0SHong Zhang ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr); 2082d6ab1dc0SHong Zhang ierr = MatSetSizes(Cmpi,pn,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 2083d6ab1dc0SHong Zhang ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr); 2084d6ab1dc0SHong Zhang ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr); 2085d6ab1dc0SHong Zhang ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 2086d6ab1dc0SHong Zhang ierr = MatSetBlockSize(Cmpi,1);CHKERRQ(ierr); 2087d6ab1dc0SHong Zhang for (i=0; i<pn; i++){ 2088d6ab1dc0SHong Zhang row = i + rstart; 2089d6ab1dc0SHong Zhang nnz = bi[i+1] - bi[i]; 2090d6ab1dc0SHong Zhang Jptr = bj + bi[i]; 2091d6ab1dc0SHong Zhang ierr = MatSetValues(Cmpi,1,&row,nnz,Jptr,vals,INSERT_VALUES);CHKERRQ(ierr); 2092d6ab1dc0SHong Zhang } 2093d6ab1dc0SHong Zhang ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2094d6ab1dc0SHong Zhang ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2095d6ab1dc0SHong Zhang ierr = PetscFree(vals);CHKERRQ(ierr); 2096d6ab1dc0SHong Zhang 2097d6ab1dc0SHong Zhang merge->bi = bi; 2098d6ab1dc0SHong Zhang merge->bj = bj; 2099d6ab1dc0SHong Zhang merge->coi = coi; 2100d6ab1dc0SHong Zhang merge->coj = coj; 2101d6ab1dc0SHong Zhang merge->buf_ri = buf_ri; 2102d6ab1dc0SHong Zhang merge->buf_rj = buf_rj; 2103d6ab1dc0SHong Zhang merge->owners_co = owners_co; 2104d6ab1dc0SHong Zhang merge->destroy = Cmpi->ops->destroy; 2105d6ab1dc0SHong Zhang merge->duplicate = Cmpi->ops->duplicate; 2106d6ab1dc0SHong Zhang 2107d6ab1dc0SHong Zhang Cmpi->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_Scalable; 2108d6ab1dc0SHong Zhang Cmpi->ops->destroy = MatDestroy_MPIAIJ_PtAP; 2109d6ab1dc0SHong Zhang 2110d6ab1dc0SHong Zhang /* attach the supporting struct to Cmpi for reuse */ 2111d6ab1dc0SHong Zhang c = (Mat_MPIAIJ*)Cmpi->data; 2112d6ab1dc0SHong Zhang c->ptap = ptap; 2113d6ab1dc0SHong Zhang ptap->api = PETSC_NULL; 2114d6ab1dc0SHong Zhang ptap->apj = PETSC_NULL; 2115d6ab1dc0SHong Zhang ptap->merge = merge; 2116d6ab1dc0SHong Zhang ptap->rmax = rmax; 2117e745005bSHong Zhang ptap->apa = PETSC_NULL; 2118f96d37fcSHong Zhang 2119f96d37fcSHong Zhang *C = Cmpi; 2120f96d37fcSHong Zhang #if defined(PETSC_USE_INFO) 2121f96d37fcSHong Zhang if (bi[pn] != 0) { 2122f96d37fcSHong Zhang ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr); 21231c7d5954SHong Zhang ierr = PetscInfo1(Cmpi,"Use MatTransposeMatMult(A,B,MatReuse,%G,&C) for best performance.\n",afill);CHKERRQ(ierr); 2124f96d37fcSHong Zhang } else { 2125f96d37fcSHong Zhang ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr); 2126f96d37fcSHong Zhang } 2127f96d37fcSHong Zhang #endif 2128f96d37fcSHong Zhang PetscFunctionReturn(0); 2129f96d37fcSHong Zhang } 2130