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> 11af0996ceSBarry Smith #include <petsc/private/vecimpl.h> 121634b032SHong Zhang 132499ec78SHong Zhang #undef __FUNCT__ 142499ec78SHong Zhang #define __FUNCT__ "MatMatMult_MPIAIJ_MPIAIJ" 152499ec78SHong Zhang PetscErrorCode MatMatMult_MPIAIJ_MPIAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill, Mat *C) 162499ec78SHong Zhang { 172499ec78SHong Zhang PetscErrorCode ierr; 180fc8cf34SHong Zhang const char *algTypes[2] = {"scalable","nonscalable"}; 190d3441aeSHong Zhang PetscInt alg=1; /* set default algorithm */ 202499ec78SHong Zhang 212499ec78SHong Zhang PetscFunctionBegin; 222499ec78SHong Zhang if (scall == MAT_INITIAL_MATRIX) { 230fc8cf34SHong Zhang ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr); 240fc8cf34SHong Zhang ierr = PetscOptionsEList("-matmatmult_via","Algorithmic approach","MatMatMult",algTypes,2,algTypes[0],&alg,NULL);CHKERRQ(ierr); 250fc8cf34SHong Zhang ierr = PetscOptionsEnd();CHKERRQ(ierr); 260fc8cf34SHong Zhang 273ff4c91cSHong Zhang ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr); 280fc8cf34SHong Zhang switch (alg) { 290fc8cf34SHong Zhang case 1: 300fc8cf34SHong Zhang ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(A,B,fill,C);CHKERRQ(ierr); 310fc8cf34SHong Zhang break; 320fc8cf34SHong Zhang default: 33b2405163SHong Zhang ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ(A,B,fill,C);CHKERRQ(ierr); 340fc8cf34SHong Zhang break; 350fc8cf34SHong Zhang } 363ff4c91cSHong Zhang ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr); 372499ec78SHong Zhang } 383ff4c91cSHong Zhang ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr); 39598bc09dSHong Zhang ierr = (*(*C)->ops->matmultnumeric)(A,B,*C);CHKERRQ(ierr); 403ff4c91cSHong Zhang ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr); 412499ec78SHong Zhang PetscFunctionReturn(0); 422499ec78SHong Zhang } 432499ec78SHong Zhang 44f33d1a9aSHong Zhang #undef __FUNCT__ 45a1a4e84aSHong Zhang #define __FUNCT__ "MatDestroy_MPIAIJ_MatMatMult" 46a1a4e84aSHong Zhang PetscErrorCode MatDestroy_MPIAIJ_MatMatMult(Mat A) 47598bc09dSHong Zhang { 48598bc09dSHong Zhang PetscErrorCode ierr; 49598bc09dSHong Zhang Mat_MPIAIJ *a = (Mat_MPIAIJ*)A->data; 50598bc09dSHong Zhang Mat_PtAPMPI *ptap = a->ptap; 51598bc09dSHong Zhang 52598bc09dSHong Zhang PetscFunctionBegin; 53b7f45c76SHong Zhang ierr = PetscFree2(ptap->startsj_s,ptap->startsj_r);CHKERRQ(ierr); 54598bc09dSHong Zhang ierr = PetscFree(ptap->bufa);CHKERRQ(ierr); 55a1a4e84aSHong Zhang ierr = MatDestroy(&ptap->P_loc);CHKERRQ(ierr); 56a1a4e84aSHong Zhang ierr = MatDestroy(&ptap->P_oth);CHKERRQ(ierr); 57ab1b013aSHong Zhang ierr = MatDestroy(&ptap->Pt);CHKERRQ(ierr); 58a1a4e84aSHong Zhang ierr = PetscFree(ptap->api);CHKERRQ(ierr); 59a1a4e84aSHong Zhang ierr = PetscFree(ptap->apj);CHKERRQ(ierr); 60598bc09dSHong Zhang ierr = PetscFree(ptap->apa);CHKERRQ(ierr); 61598bc09dSHong Zhang ierr = ptap->destroy(A);CHKERRQ(ierr); 62598bc09dSHong Zhang ierr = PetscFree(ptap);CHKERRQ(ierr); 63598bc09dSHong Zhang PetscFunctionReturn(0); 64598bc09dSHong Zhang } 65598bc09dSHong Zhang 662499ec78SHong Zhang #undef __FUNCT__ 67a1a4e84aSHong Zhang #define __FUNCT__ "MatDuplicate_MPIAIJ_MatMatMult" 68a1a4e84aSHong Zhang PetscErrorCode MatDuplicate_MPIAIJ_MatMatMult(Mat A, MatDuplicateOption op, Mat *M) 694ae0847bSHong Zhang { 704ae0847bSHong Zhang PetscErrorCode ierr; 714ae0847bSHong Zhang Mat_MPIAIJ *a = (Mat_MPIAIJ*)A->data; 724ae0847bSHong Zhang Mat_PtAPMPI *ptap = a->ptap; 734ae0847bSHong Zhang 744ae0847bSHong Zhang PetscFunctionBegin; 754ae0847bSHong Zhang ierr = (*ptap->duplicate)(A,op,M);CHKERRQ(ierr); 7626fbe8dcSKarl Rupp 774ae0847bSHong Zhang (*M)->ops->destroy = ptap->destroy; /* = MatDestroy_MPIAIJ, *M doesn't duplicate A's special structure! */ 784ae0847bSHong Zhang (*M)->ops->duplicate = ptap->duplicate; /* = MatDuplicate_MPIAIJ */ 794ae0847bSHong Zhang PetscFunctionReturn(0); 804ae0847bSHong Zhang } 814ae0847bSHong Zhang 824ae0847bSHong Zhang #undef __FUNCT__ 83f8efd71cSHong Zhang #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable" 84f8efd71cSHong Zhang PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable(Mat A,Mat P,Mat C) 85598bc09dSHong Zhang { 86598bc09dSHong Zhang PetscErrorCode ierr; 874ae0847bSHong Zhang Mat_MPIAIJ *a =(Mat_MPIAIJ*)A->data,*c=(Mat_MPIAIJ*)C->data; 88598bc09dSHong Zhang Mat_SeqAIJ *ad =(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data; 89598bc09dSHong Zhang Mat_SeqAIJ *cd =(Mat_SeqAIJ*)(c->A)->data,*co=(Mat_SeqAIJ*)(c->B)->data; 90c12557a1SHong Zhang PetscScalar *cda=cd->a,*coa=co->a; 91598bc09dSHong Zhang Mat_SeqAIJ *p_loc,*p_oth; 929ce11a7cSHong Zhang PetscScalar *apa,*ca; 93c12557a1SHong Zhang PetscInt cm =C->rmap->n; 94598bc09dSHong Zhang Mat_PtAPMPI *ptap=c->ptap; 95c12557a1SHong Zhang PetscInt *api,*apj,*apJ,i,k; 9629825010SHong Zhang PetscInt cstart=C->cmap->rstart; 97598bc09dSHong Zhang PetscInt cdnz,conz,k0,k1; 989467f45fSHong Zhang MPI_Comm comm; 999467f45fSHong Zhang PetscMPIInt size; 100598bc09dSHong Zhang 101598bc09dSHong Zhang PetscFunctionBegin; 1029467f45fSHong Zhang ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 1039467f45fSHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1049467f45fSHong Zhang 105a1a4e84aSHong Zhang /* 1) get P_oth = ptap->P_oth and P_loc = ptap->P_loc */ 106598bc09dSHong Zhang /*-----------------------------------------------------*/ 107a1a4e84aSHong Zhang /* update numerical values of P_oth and P_loc */ 108b7f45c76SHong Zhang ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 109a1a4e84aSHong Zhang ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 110598bc09dSHong Zhang 111598bc09dSHong Zhang /* 2) compute numeric C_loc = A_loc*P = Ad*P_loc + Ao*P_oth */ 112598bc09dSHong Zhang /*----------------------------------------------------------*/ 113598bc09dSHong Zhang /* get data from symbolic products */ 114a1a4e84aSHong Zhang p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 115c12557a1SHong Zhang p_oth = NULL; 1169467f45fSHong Zhang if (size >1) { 1179467f45fSHong Zhang p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 1189467f45fSHong Zhang } 119598bc09dSHong Zhang 120598bc09dSHong Zhang /* get apa for storing dense row A[i,:]*P */ 121598bc09dSHong Zhang apa = ptap->apa; 122598bc09dSHong Zhang 12329825010SHong Zhang api = ptap->api; 12429825010SHong Zhang apj = ptap->apj; 125598bc09dSHong Zhang for (i=0; i<cm; i++) { 126c12557a1SHong Zhang /* compute apa = A[i,:]*P */ 127e497d3c8SHong Zhang AProw_nonscalable(i,ad,ao,p_loc,p_oth,apa); 128598bc09dSHong Zhang 129598bc09dSHong Zhang /* set values in C */ 130598bc09dSHong Zhang apJ = apj + api[i]; 131598bc09dSHong Zhang cdnz = cd->i[i+1] - cd->i[i]; 132598bc09dSHong Zhang conz = co->i[i+1] - co->i[i]; 133598bc09dSHong Zhang 134598bc09dSHong Zhang /* 1st off-diagoanl part of C */ 135598bc09dSHong Zhang ca = coa + co->i[i]; 136598bc09dSHong Zhang k = 0; 137598bc09dSHong Zhang for (k0=0; k0<conz; k0++) { 138598bc09dSHong Zhang if (apJ[k] >= cstart) break; 139598bc09dSHong Zhang ca[k0] = apa[apJ[k]]; 140598bc09dSHong Zhang apa[apJ[k]] = 0.0; 141598bc09dSHong Zhang k++; 142598bc09dSHong Zhang } 143598bc09dSHong Zhang 144598bc09dSHong Zhang /* diagonal part of C */ 145598bc09dSHong Zhang ca = cda + cd->i[i]; 146598bc09dSHong Zhang for (k1=0; k1<cdnz; k1++) { 147598bc09dSHong Zhang ca[k1] = apa[apJ[k]]; 148598bc09dSHong Zhang apa[apJ[k]] = 0.0; 149598bc09dSHong Zhang k++; 150598bc09dSHong Zhang } 151598bc09dSHong Zhang 152598bc09dSHong Zhang /* 2nd off-diagoanl part of C */ 153598bc09dSHong Zhang ca = coa + co->i[i]; 154598bc09dSHong Zhang for (; k0<conz; k0++) { 155598bc09dSHong Zhang ca[k0] = apa[apJ[k]]; 156598bc09dSHong Zhang apa[apJ[k]] = 0.0; 157598bc09dSHong Zhang k++; 158598bc09dSHong Zhang } 159598bc09dSHong Zhang } 160598bc09dSHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 161598bc09dSHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 162598bc09dSHong Zhang PetscFunctionReturn(0); 163598bc09dSHong Zhang } 164598bc09dSHong Zhang 165598bc09dSHong Zhang #undef __FUNCT__ 1660fc8cf34SHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable" 1670fc8cf34SHong Zhang PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat A,Mat P,PetscReal fill,Mat *C) 168598bc09dSHong Zhang { 169598bc09dSHong Zhang PetscErrorCode ierr; 170ce94432eSBarry Smith MPI_Comm comm; 1719467f45fSHong Zhang PetscMPIInt size; 172bfcd1a73SHong Zhang Mat Cmpi; 173598bc09dSHong Zhang Mat_PtAPMPI *ptap; 1740298fd71SBarry Smith PetscFreeSpaceList free_space=NULL,current_space=NULL; 1754ae0847bSHong Zhang Mat_MPIAIJ *a =(Mat_MPIAIJ*)A->data,*c; 176bfcd1a73SHong Zhang Mat_SeqAIJ *ad =(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_loc,*p_oth; 1774ae0847bSHong Zhang PetscInt *pi_loc,*pj_loc,*pi_oth,*pj_oth,*dnz,*onz; 1784ae0847bSHong Zhang PetscInt *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,rstart=A->rmap->rstart; 179d14fa243SHong Zhang PetscInt *lnk,i,pnz,row,*api,*apj,*Jptr,apnz,nspacedouble=0,j,nzi; 180bfcd1a73SHong Zhang PetscInt am=A->rmap->n,pN=P->cmap->N,pn=P->cmap->n,pm=P->rmap->n; 181598bc09dSHong Zhang PetscBT lnkbt; 182598bc09dSHong Zhang PetscScalar *apa; 183bfcd1a73SHong Zhang PetscReal afill; 184f84c536bSHong Zhang PetscInt nlnk_max,armax,prmax; 185598bc09dSHong Zhang 186598bc09dSHong Zhang PetscFunctionBegin; 187ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 1889467f45fSHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1899467f45fSHong Zhang 190a1a4e84aSHong Zhang if (A->cmap->rstart != P->rmap->rstart || A->cmap->rend != P->rmap->rend) { 191cf1a0672SHong 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); 192a1a4e84aSHong Zhang } 193152983bfSHong Zhang 194598bc09dSHong Zhang /* create struct Mat_PtAPMPI and attached it to C later */ 195b00a9115SJed Brown ierr = PetscNew(&ptap);CHKERRQ(ierr); 196598bc09dSHong Zhang 197598bc09dSHong Zhang /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */ 198b7f45c76SHong Zhang ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 19919f0e57aSHong Zhang 200598bc09dSHong Zhang /* get P_loc by taking all local rows of P */ 201a1a4e84aSHong Zhang ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 202598bc09dSHong Zhang 203a1a4e84aSHong Zhang p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 204598bc09dSHong Zhang pi_loc = p_loc->i; pj_loc = p_loc->j; 2059467f45fSHong Zhang if (size > 1) { 2069467f45fSHong Zhang p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 207598bc09dSHong Zhang pi_oth = p_oth->i; pj_oth = p_oth->j; 20820e3dc75SHong Zhang } else { 20920e3dc75SHong Zhang p_oth = NULL; 21020e3dc75SHong Zhang pi_oth = NULL; pj_oth = NULL; 2119467f45fSHong Zhang } 212598bc09dSHong Zhang 213598bc09dSHong Zhang /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */ 214598bc09dSHong Zhang /*-------------------------------------------------------------------*/ 215854ce69bSBarry Smith ierr = PetscMalloc1(am+2,&api);CHKERRQ(ierr); 216a1a4e84aSHong Zhang ptap->api = api; 217598bc09dSHong Zhang api[0] = 0; 218598bc09dSHong Zhang 219598bc09dSHong Zhang /* create and initialize a linked list */ 220f84c536bSHong Zhang armax = ad->rmax+ao->rmax; 2219467f45fSHong Zhang if (size >1) { 222f84c536bSHong Zhang prmax = PetscMax(p_loc->rmax,p_oth->rmax); 2239467f45fSHong Zhang } else { 2249467f45fSHong Zhang prmax = p_loc->rmax; 2259467f45fSHong Zhang } 226*e5541957SHong Zhang nlnk_max = armax*prmax; /* not good -- needs fix!!! */ 227f84c536bSHong Zhang if (!nlnk_max || nlnk_max > pN) nlnk_max = pN; 2280d3134e9SHong Zhang ierr = PetscLLCondensedCreate(nlnk_max,pN,&lnk,&lnkbt);CHKERRQ(ierr); 229598bc09dSHong Zhang 230bfcd1a73SHong Zhang /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */ 231bfcd1a73SHong Zhang ierr = PetscFreeSpaceGet((PetscInt)(fill*(adi[am]+aoi[am]+pi_loc[pm])),&free_space);CHKERRQ(ierr); 2322205254eSKarl Rupp 233598bc09dSHong Zhang current_space = free_space; 234598bc09dSHong Zhang 235598bc09dSHong Zhang ierr = MatPreallocateInitialize(comm,am,pn,dnz,onz);CHKERRQ(ierr); 236598bc09dSHong Zhang for (i=0; i<am; i++) { 237598bc09dSHong Zhang /* diagonal portion of A */ 238598bc09dSHong Zhang nzi = adi[i+1] - adi[i]; 239598bc09dSHong Zhang for (j=0; j<nzi; j++) { 240598bc09dSHong Zhang row = *adj++; 241598bc09dSHong Zhang pnz = pi_loc[row+1] - pi_loc[row]; 242598bc09dSHong Zhang Jptr = pj_loc + pi_loc[row]; 243598bc09dSHong Zhang /* add non-zero cols of P into the sorted linked list lnk */ 2441fbbb359SHong Zhang ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 245598bc09dSHong Zhang } 246598bc09dSHong Zhang /* off-diagonal portion of A */ 247598bc09dSHong Zhang nzi = aoi[i+1] - aoi[i]; 248598bc09dSHong Zhang for (j=0; j<nzi; j++) { 249598bc09dSHong Zhang row = *aoj++; 250598bc09dSHong Zhang pnz = pi_oth[row+1] - pi_oth[row]; 251598bc09dSHong Zhang Jptr = pj_oth + pi_oth[row]; 2521fbbb359SHong Zhang ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 253598bc09dSHong Zhang } 254598bc09dSHong Zhang 255d14fa243SHong Zhang apnz = lnk[0]; 256598bc09dSHong Zhang api[i+1] = api[i] + apnz; 257598bc09dSHong Zhang 258598bc09dSHong Zhang /* if free space is not available, double the total space in the list */ 259598bc09dSHong Zhang if (current_space->local_remaining<apnz) { 260598bc09dSHong Zhang ierr = PetscFreeSpaceGet(apnz+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 261598bc09dSHong Zhang nspacedouble++; 262598bc09dSHong Zhang } 263598bc09dSHong Zhang 264598bc09dSHong Zhang /* Copy data into free space, then initialize lnk */ 2651fbbb359SHong Zhang ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 266598bc09dSHong Zhang ierr = MatPreallocateSet(i+rstart,apnz,current_space->array,dnz,onz);CHKERRQ(ierr); 2672205254eSKarl Rupp 268598bc09dSHong Zhang current_space->array += apnz; 269598bc09dSHong Zhang current_space->local_used += apnz; 270598bc09dSHong Zhang current_space->local_remaining -= apnz; 271598bc09dSHong Zhang } 272598bc09dSHong Zhang 273598bc09dSHong Zhang /* Allocate space for apj, initialize apj, and */ 274598bc09dSHong Zhang /* destroy list of free space and other temporary array(s) */ 275854ce69bSBarry Smith ierr = PetscMalloc1(api[am]+1,&ptap->apj);CHKERRQ(ierr); 276a1a4e84aSHong Zhang apj = ptap->apj; 277a1a4e84aSHong Zhang ierr = PetscFreeSpaceContiguous(&free_space,ptap->apj);CHKERRQ(ierr); 278598bc09dSHong Zhang ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 279598bc09dSHong Zhang 280f84c536bSHong Zhang /* malloc apa to store dense row A[i,:]*P */ 2811795a4d1SJed Brown ierr = PetscCalloc1(pN,&apa);CHKERRQ(ierr); 2822205254eSKarl Rupp 283f84c536bSHong Zhang ptap->apa = apa; 284f84c536bSHong Zhang 285bfcd1a73SHong Zhang /* create and assemble symbolic parallel matrix Cmpi */ 286598bc09dSHong Zhang /*----------------------------------------------------*/ 287bfcd1a73SHong Zhang ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr); 288bfcd1a73SHong Zhang ierr = MatSetSizes(Cmpi,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 28933d57670SJed Brown ierr = MatSetBlockSizesFromMats(Cmpi,A,P);CHKERRQ(ierr); 290a2f3521dSMark F. Adams 291bfcd1a73SHong Zhang ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr); 292bfcd1a73SHong Zhang ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr); 293598bc09dSHong Zhang ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 294598bc09dSHong Zhang for (i=0; i<am; i++) { 295598bc09dSHong Zhang row = i + rstart; 296598bc09dSHong Zhang apnz = api[i+1] - api[i]; 297bfcd1a73SHong Zhang ierr = MatSetValues(Cmpi,1,&row,apnz,apj,apa,INSERT_VALUES);CHKERRQ(ierr); 298598bc09dSHong Zhang apj += apnz; 299598bc09dSHong Zhang } 300bfcd1a73SHong Zhang ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 301bfcd1a73SHong Zhang ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 302598bc09dSHong Zhang 303bfcd1a73SHong Zhang ptap->destroy = Cmpi->ops->destroy; 304bfcd1a73SHong Zhang ptap->duplicate = Cmpi->ops->duplicate; 305f8efd71cSHong Zhang Cmpi->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable; 306bfcd1a73SHong Zhang Cmpi->ops->destroy = MatDestroy_MPIAIJ_MatMatMult; 307bfcd1a73SHong Zhang Cmpi->ops->duplicate = MatDuplicate_MPIAIJ_MatMatMult; 308598bc09dSHong Zhang 309bfcd1a73SHong Zhang /* attach the supporting struct to Cmpi for reuse */ 310bfcd1a73SHong Zhang c = (Mat_MPIAIJ*)Cmpi->data; 311598bc09dSHong Zhang c->ptap = ptap; 312598bc09dSHong Zhang 313bfcd1a73SHong Zhang *C = Cmpi; 314bfcd1a73SHong Zhang 315bfcd1a73SHong Zhang /* set MatInfo */ 316118727c9SMark F. Adams afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1) + 1.e-5; 317bfcd1a73SHong Zhang if (afill < 1.0) afill = 1.0; 318bfcd1a73SHong Zhang Cmpi->info.mallocs = nspacedouble; 319bfcd1a73SHong Zhang Cmpi->info.fill_ratio_given = fill; 320bfcd1a73SHong Zhang Cmpi->info.fill_ratio_needed = afill; 321bfcd1a73SHong Zhang 322bfcd1a73SHong Zhang #if defined(PETSC_USE_INFO) 323bfcd1a73SHong Zhang if (api[am]) { 32457622a8eSBarry Smith ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr); 32557622a8eSBarry Smith ierr = PetscInfo1(Cmpi,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr); 326bfcd1a73SHong Zhang } else { 327bfcd1a73SHong Zhang ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr); 328bfcd1a73SHong Zhang } 329bfcd1a73SHong Zhang #endif 330598bc09dSHong Zhang PetscFunctionReturn(0); 331598bc09dSHong Zhang } 332598bc09dSHong Zhang 3339123193aSHong Zhang #undef __FUNCT__ 3349123193aSHong Zhang #define __FUNCT__ "MatMatMult_MPIAIJ_MPIDense" 3359123193aSHong Zhang PetscErrorCode MatMatMult_MPIAIJ_MPIDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 3369123193aSHong Zhang { 3379123193aSHong Zhang PetscErrorCode ierr; 3389123193aSHong Zhang 3399123193aSHong Zhang PetscFunctionBegin; 3409123193aSHong Zhang if (scall == MAT_INITIAL_MATRIX) { 3413ff4c91cSHong Zhang ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr); 3429123193aSHong Zhang ierr = MatMatMultSymbolic_MPIAIJ_MPIDense(A,B,fill,C);CHKERRQ(ierr); 3433ff4c91cSHong Zhang ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr); 3449123193aSHong Zhang } 3453ff4c91cSHong Zhang ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr); 3469123193aSHong Zhang ierr = MatMatMultNumeric_MPIAIJ_MPIDense(A,B,*C);CHKERRQ(ierr); 3473ff4c91cSHong Zhang ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr); 3489123193aSHong Zhang PetscFunctionReturn(0); 3499123193aSHong Zhang } 3509123193aSHong Zhang 3514324174eSBarry Smith typedef struct { 3524324174eSBarry Smith Mat workB; 3534324174eSBarry Smith PetscScalar *rvalues,*svalues; 3544324174eSBarry Smith MPI_Request *rwaits,*swaits; 3554324174eSBarry Smith } MPIAIJ_MPIDense; 3564324174eSBarry Smith 3577af9e4fdSHong Zhang #undef __FUNCT__ 35896b95a6bSBarry Smith #define __FUNCT__ "MatMPIAIJ_MPIDenseDestroy" 35996b95a6bSBarry Smith PetscErrorCode MatMPIAIJ_MPIDenseDestroy(void *ctx) 3604324174eSBarry Smith { 3614324174eSBarry Smith MPIAIJ_MPIDense *contents = (MPIAIJ_MPIDense*) ctx; 3624324174eSBarry Smith PetscErrorCode ierr; 3634324174eSBarry Smith 3644324174eSBarry Smith PetscFunctionBegin; 3656bf464f9SBarry Smith ierr = MatDestroy(&contents->workB);CHKERRQ(ierr); 3664324174eSBarry Smith ierr = PetscFree4(contents->rvalues,contents->svalues,contents->rwaits,contents->swaits);CHKERRQ(ierr); 3674324174eSBarry Smith ierr = PetscFree(contents);CHKERRQ(ierr); 3684324174eSBarry Smith PetscFunctionReturn(0); 3694324174eSBarry Smith } 3704324174eSBarry Smith 3719123193aSHong Zhang #undef __FUNCT__ 372fd4e9aacSBarry Smith #define __FUNCT__ "MatMatMultNumeric_MPIDense" 373fd4e9aacSBarry Smith /* 374fd4e9aacSBarry Smith This is a "dummy function" that handles the case where matrix C was created as a dense matrix 375fd4e9aacSBarry Smith directly by the user and passed to MatMatMult() with the MAT_REUSE_MATRIX option 376fd4e9aacSBarry Smith 377fd4e9aacSBarry Smith It is the same as MatMatMultSymbolic_MPIAIJ_MPIDense() except does not create C 378fd4e9aacSBarry Smith */ 379fd4e9aacSBarry Smith PetscErrorCode MatMatMultNumeric_MPIDense(Mat A,Mat B,Mat C) 380fd4e9aacSBarry Smith { 381fd4e9aacSBarry Smith PetscErrorCode ierr; 382fd4e9aacSBarry Smith PetscBool flg; 383fd4e9aacSBarry Smith Mat_MPIAIJ *aij = (Mat_MPIAIJ*) A->data; 384fd4e9aacSBarry Smith PetscInt nz = aij->B->cmap->n; 385fd4e9aacSBarry Smith PetscContainer container; 386fd4e9aacSBarry Smith MPIAIJ_MPIDense *contents; 387fd4e9aacSBarry Smith VecScatter ctx = aij->Mvctx; 388fd4e9aacSBarry Smith VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata; 389fd4e9aacSBarry Smith VecScatter_MPI_General *to = (VecScatter_MPI_General*) ctx->todata; 390fd4e9aacSBarry Smith 391fd4e9aacSBarry Smith PetscFunctionBegin; 392fd4e9aacSBarry Smith ierr = PetscObjectTypeCompare((PetscObject)B,MATMPIDENSE,&flg);CHKERRQ(ierr); 393fd4e9aacSBarry Smith if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Second matrix must be mpidense"); 394fd4e9aacSBarry Smith 395fd4e9aacSBarry Smith /* Handle case where where user provided the final C matrix rather than calling MatMatMult() with MAT_INITIAL_MATRIX*/ 396fd4e9aacSBarry Smith ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&flg);CHKERRQ(ierr); 397fd4e9aacSBarry Smith if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"First matrix must be MPIAIJ"); 398fd4e9aacSBarry Smith 399fd4e9aacSBarry Smith C->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIDense; 400fd4e9aacSBarry Smith 401fd4e9aacSBarry Smith ierr = PetscNew(&contents);CHKERRQ(ierr); 402fd4e9aacSBarry Smith /* Create work matrix used to store off processor rows of B needed for local product */ 403fd4e9aacSBarry Smith ierr = MatCreateSeqDense(PETSC_COMM_SELF,nz,B->cmap->N,NULL,&contents->workB);CHKERRQ(ierr); 404fd4e9aacSBarry Smith /* Create work arrays needed */ 405fd4e9aacSBarry Smith ierr = PetscMalloc4(B->cmap->N*from->starts[from->n],&contents->rvalues, 406fd4e9aacSBarry Smith B->cmap->N*to->starts[to->n],&contents->svalues, 407fd4e9aacSBarry Smith from->n,&contents->rwaits, 408fd4e9aacSBarry Smith to->n,&contents->swaits);CHKERRQ(ierr); 409fd4e9aacSBarry Smith 410fd4e9aacSBarry Smith ierr = PetscContainerCreate(PetscObjectComm((PetscObject)A),&container);CHKERRQ(ierr); 411fd4e9aacSBarry Smith ierr = PetscContainerSetPointer(container,contents);CHKERRQ(ierr); 412fd4e9aacSBarry Smith ierr = PetscContainerSetUserDestroy(container,MatMPIAIJ_MPIDenseDestroy);CHKERRQ(ierr); 413fd4e9aacSBarry Smith ierr = PetscObjectCompose((PetscObject)C,"workB",(PetscObject)container);CHKERRQ(ierr); 414fd4e9aacSBarry Smith ierr = PetscContainerDestroy(&container);CHKERRQ(ierr); 415fd4e9aacSBarry Smith 416fd4e9aacSBarry Smith ierr = (*C->ops->matmultnumeric)(A,B,C);CHKERRQ(ierr); 417fd4e9aacSBarry Smith PetscFunctionReturn(0); 418fd4e9aacSBarry Smith } 419fd4e9aacSBarry Smith 420fd4e9aacSBarry Smith #undef __FUNCT__ 4219123193aSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIDense" 4229123193aSHong Zhang PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIDense(Mat A,Mat B,PetscReal fill,Mat *C) 4239123193aSHong Zhang { 4249123193aSHong Zhang PetscErrorCode ierr; 425f32d5d43SBarry Smith Mat_MPIAIJ *aij = (Mat_MPIAIJ*) A->data; 426d0f46423SBarry Smith PetscInt nz = aij->B->cmap->n; 427bf0cc555SLisandro Dalcin PetscContainer container; 4284324174eSBarry Smith MPIAIJ_MPIDense *contents; 4294324174eSBarry Smith VecScatter ctx = aij->Mvctx; 4304324174eSBarry Smith VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata; 4314324174eSBarry Smith VecScatter_MPI_General *to = (VecScatter_MPI_General*) ctx->todata; 432d0f46423SBarry Smith PetscInt m = A->rmap->n,n=B->cmap->n; 4339123193aSHong Zhang 4349123193aSHong Zhang PetscFunctionBegin; 435ce94432eSBarry Smith ierr = MatCreate(PetscObjectComm((PetscObject)B),C);CHKERRQ(ierr); 436d0f46423SBarry Smith ierr = MatSetSizes(*C,m,n,A->rmap->N,B->cmap->N);CHKERRQ(ierr); 43733d57670SJed Brown ierr = MatSetBlockSizesFromMats(*C,A,B);CHKERRQ(ierr); 438cb2480eaSBarry Smith ierr = MatSetType(*C,MATMPIDENSE);CHKERRQ(ierr); 4390298fd71SBarry Smith ierr = MatMPIDenseSetPreallocation(*C,NULL);CHKERRQ(ierr); 440cb2480eaSBarry Smith ierr = MatAssemblyBegin(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 441cb2480eaSBarry Smith ierr = MatAssemblyEnd(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 4422205254eSKarl Rupp 443d73949e8SHong Zhang (*C)->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIDense; 444f32d5d43SBarry Smith 445b00a9115SJed Brown ierr = PetscNew(&contents);CHKERRQ(ierr); 446f32d5d43SBarry Smith /* Create work matrix used to store off processor rows of B needed for local product */ 4470298fd71SBarry Smith ierr = MatCreateSeqDense(PETSC_COMM_SELF,nz,B->cmap->N,NULL,&contents->workB);CHKERRQ(ierr); 4484324174eSBarry Smith /* Create work arrays needed */ 449dcca6d9dSJed Brown ierr = PetscMalloc4(B->cmap->N*from->starts[from->n],&contents->rvalues, 450dcca6d9dSJed Brown B->cmap->N*to->starts[to->n],&contents->svalues, 451dcca6d9dSJed Brown from->n,&contents->rwaits, 452dcca6d9dSJed Brown to->n,&contents->swaits);CHKERRQ(ierr); 4534324174eSBarry Smith 454ce94432eSBarry Smith ierr = PetscContainerCreate(PetscObjectComm((PetscObject)A),&container);CHKERRQ(ierr); 455bf0cc555SLisandro Dalcin ierr = PetscContainerSetPointer(container,contents);CHKERRQ(ierr); 45696b95a6bSBarry Smith ierr = PetscContainerSetUserDestroy(container,MatMPIAIJ_MPIDenseDestroy);CHKERRQ(ierr); 457bf0cc555SLisandro Dalcin ierr = PetscObjectCompose((PetscObject)(*C),"workB",(PetscObject)container);CHKERRQ(ierr); 458bf0cc555SLisandro Dalcin ierr = PetscContainerDestroy(&container);CHKERRQ(ierr); 4599123193aSHong Zhang PetscFunctionReturn(0); 4609123193aSHong Zhang } 4619123193aSHong Zhang 4627af9e4fdSHong Zhang #undef __FUNCT__ 4637af9e4fdSHong Zhang #define __FUNCT__ "MatMPIDenseScatter" 464f32d5d43SBarry Smith /* 4652636ff26SBarry Smith Performs an efficient scatter on the rows of B needed by this process; this is 4662636ff26SBarry Smith a modification of the VecScatterBegin_() routines. 467f32d5d43SBarry Smith */ 4684324174eSBarry Smith PetscErrorCode MatMPIDenseScatter(Mat A,Mat B,Mat C,Mat *outworkB) 469f32d5d43SBarry Smith { 470f32d5d43SBarry Smith Mat_MPIAIJ *aij = (Mat_MPIAIJ*)A->data; 471f32d5d43SBarry Smith PetscErrorCode ierr; 472f32d5d43SBarry Smith PetscScalar *b,*w,*svalues,*rvalues; 473f32d5d43SBarry Smith VecScatter ctx = aij->Mvctx; 474f32d5d43SBarry Smith VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata; 475f32d5d43SBarry Smith VecScatter_MPI_General *to = (VecScatter_MPI_General*) ctx->todata; 476f32d5d43SBarry Smith PetscInt i,j,k; 477aa5bb8c0SSatish Balay PetscInt *sindices,*sstarts,*rindices,*rstarts; 478aa5bb8c0SSatish Balay PetscMPIInt *sprocs,*rprocs,nrecvs; 479f32d5d43SBarry Smith MPI_Request *swaits,*rwaits; 480ce94432eSBarry Smith MPI_Comm comm; 481d0f46423SBarry Smith PetscMPIInt tag = ((PetscObject)ctx)->tag,ncols = B->cmap->N, nrows = aij->B->cmap->n,imdex,nrowsB = B->rmap->n; 482f32d5d43SBarry Smith MPI_Status status; 4834324174eSBarry Smith MPIAIJ_MPIDense *contents; 484bf0cc555SLisandro Dalcin PetscContainer container; 4854324174eSBarry Smith Mat workB; 486f32d5d43SBarry Smith 487f32d5d43SBarry Smith PetscFunctionBegin; 488ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 489bf0cc555SLisandro Dalcin ierr = PetscObjectQuery((PetscObject)C,"workB",(PetscObject*)&container);CHKERRQ(ierr); 49029825010SHong Zhang if (!container) SETERRQ(comm,PETSC_ERR_PLIB,"Container does not exist"); 491bf0cc555SLisandro Dalcin ierr = PetscContainerGetPointer(container,(void**)&contents);CHKERRQ(ierr); 4924324174eSBarry Smith 4934324174eSBarry Smith workB = *outworkB = contents->workB; 494cf1a0672SHong 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); 495f32d5d43SBarry Smith sindices = to->indices; 496f32d5d43SBarry Smith sstarts = to->starts; 497f32d5d43SBarry Smith sprocs = to->procs; 4984324174eSBarry Smith swaits = contents->swaits; 4994324174eSBarry Smith svalues = contents->svalues; 500f32d5d43SBarry Smith 501f32d5d43SBarry Smith rindices = from->indices; 502f32d5d43SBarry Smith rstarts = from->starts; 503f32d5d43SBarry Smith rprocs = from->procs; 5044324174eSBarry Smith rwaits = contents->rwaits; 5054324174eSBarry Smith rvalues = contents->rvalues; 506f32d5d43SBarry Smith 5078c778c55SBarry Smith ierr = MatDenseGetArray(B,&b);CHKERRQ(ierr); 5088c778c55SBarry Smith ierr = MatDenseGetArray(workB,&w);CHKERRQ(ierr); 509f32d5d43SBarry Smith 510f32d5d43SBarry Smith for (i=0; i<from->n; i++) { 511f32d5d43SBarry Smith ierr = MPI_Irecv(rvalues+ncols*rstarts[i],ncols*(rstarts[i+1]-rstarts[i]),MPIU_SCALAR,rprocs[i],tag,comm,rwaits+i);CHKERRQ(ierr); 512f32d5d43SBarry Smith } 513f32d5d43SBarry Smith 514f32d5d43SBarry Smith for (i=0; i<to->n; i++) { 515f32d5d43SBarry Smith /* pack a message at a time */ 516f32d5d43SBarry Smith for (j=0; j<sstarts[i+1]-sstarts[i]; j++) { 517f32d5d43SBarry Smith for (k=0; k<ncols; k++) { 518f32d5d43SBarry Smith svalues[ncols*(sstarts[i] + j) + k] = b[sindices[sstarts[i]+j] + nrowsB*k]; 5192636ff26SBarry Smith } 5202636ff26SBarry Smith } 521f32d5d43SBarry Smith ierr = MPI_Isend(svalues+ncols*sstarts[i],ncols*(sstarts[i+1]-sstarts[i]),MPIU_SCALAR,sprocs[i],tag,comm,swaits+i);CHKERRQ(ierr); 522f32d5d43SBarry Smith } 5232636ff26SBarry Smith 524f32d5d43SBarry Smith nrecvs = from->n; 525f32d5d43SBarry Smith while (nrecvs) { 526f32d5d43SBarry Smith ierr = MPI_Waitany(from->n,rwaits,&imdex,&status);CHKERRQ(ierr); 527f32d5d43SBarry Smith nrecvs--; 528f32d5d43SBarry Smith /* unpack a message at a time */ 529f32d5d43SBarry Smith for (j=0; j<rstarts[imdex+1]-rstarts[imdex]; j++) { 530f32d5d43SBarry Smith for (k=0; k<ncols; k++) { 531f32d5d43SBarry Smith w[rindices[rstarts[imdex]+j] + nrows*k] = rvalues[ncols*(rstarts[imdex] + j) + k]; 5322636ff26SBarry Smith } 5332636ff26SBarry Smith } 534f32d5d43SBarry Smith } 535cb9801acSJed Brown if (to->n) {ierr = MPI_Waitall(to->n,swaits,to->sstatus);CHKERRQ(ierr);} 536f32d5d43SBarry Smith 5378c778c55SBarry Smith ierr = MatDenseRestoreArray(B,&b);CHKERRQ(ierr); 5388c778c55SBarry Smith ierr = MatDenseRestoreArray(workB,&w);CHKERRQ(ierr); 539f32d5d43SBarry Smith ierr = MatAssemblyBegin(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 540f32d5d43SBarry Smith ierr = MatAssemblyEnd(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 541f32d5d43SBarry Smith PetscFunctionReturn(0); 542f32d5d43SBarry Smith } 543f32d5d43SBarry Smith extern PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat,Mat,Mat); 544f32d5d43SBarry Smith 5459123193aSHong Zhang #undef __FUNCT__ 5469123193aSHong Zhang #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIDense" 5479123193aSHong Zhang PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIDense(Mat A,Mat B,Mat C) 5489123193aSHong Zhang { 5499123193aSHong Zhang PetscErrorCode ierr; 550f32d5d43SBarry Smith Mat_MPIAIJ *aij = (Mat_MPIAIJ*)A->data; 551f32d5d43SBarry Smith Mat_MPIDense *bdense = (Mat_MPIDense*)B->data; 552f32d5d43SBarry Smith Mat_MPIDense *cdense = (Mat_MPIDense*)C->data; 553f32d5d43SBarry Smith Mat workB; 5549123193aSHong Zhang 5559123193aSHong Zhang PetscFunctionBegin; 556f32d5d43SBarry Smith /* diagonal block of A times all local rows of B*/ 557f32d5d43SBarry Smith ierr = MatMatMultNumeric_SeqAIJ_SeqDense(aij->A,bdense->A,cdense->A);CHKERRQ(ierr); 558f32d5d43SBarry Smith 559f32d5d43SBarry Smith /* get off processor parts of B needed to complete the product */ 5604324174eSBarry Smith ierr = MatMPIDenseScatter(A,B,C,&workB);CHKERRQ(ierr); 561f32d5d43SBarry Smith 562f32d5d43SBarry Smith /* off-diagonal block of A times nonlocal rows of B */ 563f32d5d43SBarry Smith ierr = MatMatMultNumericAdd_SeqAIJ_SeqDense(aij->B,workB,cdense->A);CHKERRQ(ierr); 5649123193aSHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 5659123193aSHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 5669123193aSHong Zhang PetscFunctionReturn(0); 5679123193aSHong Zhang } 568cf1a0672SHong Zhang 5691634b032SHong Zhang #undef __FUNCT__ 570f8efd71cSHong Zhang #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIAIJ" 571f8efd71cSHong Zhang PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C) 5721634b032SHong Zhang { 5731634b032SHong Zhang PetscErrorCode ierr; 574cf1a0672SHong Zhang Mat_MPIAIJ *a = (Mat_MPIAIJ*)A->data,*c=(Mat_MPIAIJ*)C->data; 575cf1a0672SHong Zhang Mat_SeqAIJ *ad = (Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data; 576cf1a0672SHong Zhang Mat_SeqAIJ *cd = (Mat_SeqAIJ*)(c->A)->data,*co=(Mat_SeqAIJ*)(c->B)->data; 577cf1a0672SHong Zhang PetscInt *adi = ad->i,*adj,*aoi=ao->i,*aoj; 578cf1a0672SHong Zhang PetscScalar *ada,*aoa,*cda=cd->a,*coa=co->a; 579cf1a0672SHong Zhang Mat_SeqAIJ *p_loc,*p_oth; 580cf1a0672SHong Zhang PetscInt *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pj; 58129825010SHong Zhang PetscScalar *pa_loc,*pa_oth,*pa,valtmp,*ca; 58229825010SHong Zhang PetscInt cm = C->rmap->n,anz,pnz; 583cf1a0672SHong Zhang Mat_PtAPMPI *ptap = c->ptap; 58429825010SHong Zhang PetscScalar *apa_sparse = ptap->apa; 585cf1a0672SHong Zhang PetscInt *api,*apj,*apJ,i,j,k,row; 58629825010SHong Zhang PetscInt cstart = C->cmap->rstart; 58729825010SHong Zhang PetscInt cdnz,conz,k0,k1,nextp; 588a7c7454dSHong Zhang MPI_Comm comm; 589a7c7454dSHong Zhang PetscMPIInt size; 5901634b032SHong Zhang 5911634b032SHong Zhang PetscFunctionBegin; 592a7c7454dSHong Zhang ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 593a7c7454dSHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 594a7c7454dSHong Zhang 595cf1a0672SHong Zhang /* 1) get P_oth = ptap->P_oth and P_loc = ptap->P_loc */ 596cf1a0672SHong Zhang /*-----------------------------------------------------*/ 597cf1a0672SHong Zhang /* update numerical values of P_oth and P_loc */ 598b7f45c76SHong Zhang ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 599cf1a0672SHong Zhang ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 600cf1a0672SHong Zhang 601cf1a0672SHong Zhang /* 2) compute numeric C_loc = A_loc*P = Ad*P_loc + Ao*P_oth */ 602cf1a0672SHong Zhang /*----------------------------------------------------------*/ 603cf1a0672SHong Zhang /* get data from symbolic products */ 604cf1a0672SHong Zhang p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 605cf1a0672SHong Zhang pi_loc = p_loc->i; pj_loc = p_loc->j; pa_loc = p_loc->a; 606a7c7454dSHong Zhang if (size >1) { 607a7c7454dSHong Zhang p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 608cf1a0672SHong Zhang pi_oth = p_oth->i; pj_oth = p_oth->j; pa_oth = p_oth->a; 60920e3dc75SHong Zhang } else { 61020e3dc75SHong Zhang p_oth = NULL; pi_oth = NULL; pj_oth = NULL; pa_oth = NULL; 611a7c7454dSHong Zhang } 612cf1a0672SHong Zhang 61329825010SHong Zhang api = ptap->api; 61429825010SHong Zhang apj = ptap->apj; 615cf1a0672SHong Zhang for (i=0; i<cm; i++) { 61629825010SHong Zhang apJ = apj + api[i]; 61729825010SHong Zhang 618cf1a0672SHong Zhang /* diagonal portion of A */ 619cf1a0672SHong Zhang anz = adi[i+1] - adi[i]; 620cf1a0672SHong Zhang adj = ad->j + adi[i]; 621cf1a0672SHong Zhang ada = ad->a + adi[i]; 622cf1a0672SHong Zhang for (j=0; j<anz; j++) { 623cf1a0672SHong Zhang row = adj[j]; 624cf1a0672SHong Zhang pnz = pi_loc[row+1] - pi_loc[row]; 625cf1a0672SHong Zhang pj = pj_loc + pi_loc[row]; 626cf1a0672SHong Zhang pa = pa_loc + pi_loc[row]; 62729825010SHong Zhang /* perform sparse axpy */ 628cf1a0672SHong Zhang valtmp = ada[j]; 62929825010SHong Zhang nextp = 0; 63029825010SHong Zhang for (k=0; nextp<pnz; k++) { 63129825010SHong Zhang if (apJ[k] == pj[nextp]) { /* column of AP == column of P */ 63229825010SHong Zhang apa_sparse[k] += valtmp*pa[nextp++]; 63329825010SHong Zhang } 6341634b032SHong Zhang } 635cf1a0672SHong Zhang ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr); 636cf1a0672SHong Zhang } 6371634b032SHong Zhang 638cf1a0672SHong Zhang /* off-diagonal portion of A */ 639cf1a0672SHong Zhang anz = aoi[i+1] - aoi[i]; 640cf1a0672SHong Zhang aoj = ao->j + aoi[i]; 641cf1a0672SHong Zhang aoa = ao->a + aoi[i]; 642cf1a0672SHong Zhang for (j=0; j<anz; j++) { 643cf1a0672SHong Zhang row = aoj[j]; 644cf1a0672SHong Zhang pnz = pi_oth[row+1] - pi_oth[row]; 645cf1a0672SHong Zhang pj = pj_oth + pi_oth[row]; 646cf1a0672SHong Zhang pa = pa_oth + pi_oth[row]; 64729825010SHong Zhang /* perform sparse axpy */ 648cf1a0672SHong Zhang valtmp = aoa[j]; 64929825010SHong Zhang nextp = 0; 65029825010SHong Zhang for (k=0; nextp<pnz; k++) { 65129825010SHong Zhang if (apJ[k] == pj[nextp]) { /* column of AP == column of P */ 65229825010SHong Zhang apa_sparse[k] += valtmp*pa[nextp++]; 65329825010SHong Zhang } 654cf1a0672SHong Zhang } 655cf1a0672SHong Zhang ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr); 656cf1a0672SHong Zhang } 6571634b032SHong Zhang 658cf1a0672SHong Zhang /* set values in C */ 659cf1a0672SHong Zhang cdnz = cd->i[i+1] - cd->i[i]; 660cf1a0672SHong Zhang conz = co->i[i+1] - co->i[i]; 6611634b032SHong Zhang 662cf1a0672SHong Zhang /* 1st off-diagoanl part of C */ 663cf1a0672SHong Zhang ca = coa + co->i[i]; 664cf1a0672SHong Zhang k = 0; 665cf1a0672SHong Zhang for (k0=0; k0<conz; k0++) { 666cf1a0672SHong Zhang if (apJ[k] >= cstart) break; 66729825010SHong Zhang ca[k0] = apa_sparse[k]; 66829825010SHong Zhang apa_sparse[k] = 0.0; 669cf1a0672SHong Zhang k++; 670cf1a0672SHong Zhang } 6711634b032SHong Zhang 672cf1a0672SHong Zhang /* diagonal part of C */ 673cf1a0672SHong Zhang ca = cda + cd->i[i]; 674cf1a0672SHong Zhang for (k1=0; k1<cdnz; k1++) { 67529825010SHong Zhang ca[k1] = apa_sparse[k]; 67629825010SHong Zhang apa_sparse[k] = 0.0; 677cf1a0672SHong Zhang k++; 678cf1a0672SHong Zhang } 679cf1a0672SHong Zhang 680cf1a0672SHong Zhang /* 2nd off-diagoanl part of C */ 681cf1a0672SHong Zhang ca = coa + co->i[i]; 682cf1a0672SHong Zhang for (; k0<conz; k0++) { 68329825010SHong Zhang ca[k0] = apa_sparse[k]; 68429825010SHong Zhang apa_sparse[k] = 0.0; 685cf1a0672SHong Zhang k++; 686cf1a0672SHong Zhang } 687cf1a0672SHong Zhang } 688cf1a0672SHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 689cf1a0672SHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 690cf1a0672SHong Zhang PetscFunctionReturn(0); 691cf1a0672SHong Zhang } 692cf1a0672SHong Zhang 6930fc8cf34SHong Zhang /* same as MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(), except using LLCondensed to avoid O(BN) memory requirement */ 694cf1a0672SHong Zhang #undef __FUNCT__ 695b2405163SHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ" 696b2405163SHong Zhang PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C) 697cf1a0672SHong Zhang { 698cf1a0672SHong Zhang PetscErrorCode ierr; 699ce94432eSBarry Smith MPI_Comm comm; 700a7c7454dSHong Zhang PetscMPIInt size; 701cf1a0672SHong Zhang Mat Cmpi; 702cf1a0672SHong Zhang Mat_PtAPMPI *ptap; 7030298fd71SBarry Smith PetscFreeSpaceList free_space = NULL,current_space=NULL; 704cf1a0672SHong Zhang Mat_MPIAIJ *a = (Mat_MPIAIJ*)A->data,*c; 705cf1a0672SHong Zhang Mat_SeqAIJ *ad = (Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_loc,*p_oth; 706cf1a0672SHong Zhang PetscInt *pi_loc,*pj_loc,*pi_oth,*pj_oth,*dnz,*onz; 707cf1a0672SHong Zhang PetscInt *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,rstart=A->rmap->rstart; 708f84c536bSHong Zhang PetscInt i,pnz,row,*api,*apj,*Jptr,apnz,nspacedouble=0,j,nzi,*lnk,apnz_max=0; 709cf1a0672SHong Zhang PetscInt am=A->rmap->n,pN=P->cmap->N,pn=P->cmap->n,pm=P->rmap->n; 710cf1a0672SHong Zhang PetscInt nlnk_max,armax,prmax; 711cf1a0672SHong Zhang PetscReal afill; 712cf1a0672SHong Zhang PetscScalar *apa; 713cf1a0672SHong Zhang 714cf1a0672SHong Zhang PetscFunctionBegin; 715ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 716a7c7454dSHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 717a7c7454dSHong Zhang 718cf1a0672SHong Zhang /* create struct Mat_PtAPMPI and attached it to C later */ 719b00a9115SJed Brown ierr = PetscNew(&ptap);CHKERRQ(ierr); 720cf1a0672SHong Zhang 721cf1a0672SHong Zhang /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */ 722b7f45c76SHong Zhang ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 72319f0e57aSHong Zhang 724cf1a0672SHong Zhang /* get P_loc by taking all local rows of P */ 725cf1a0672SHong Zhang ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 726cf1a0672SHong Zhang 727cf1a0672SHong Zhang p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 728cf1a0672SHong Zhang pi_loc = p_loc->i; pj_loc = p_loc->j; 729a7c7454dSHong Zhang if (size > 1) { 730a7c7454dSHong Zhang p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 73120e3dc75SHong Zhang pi_oth = p_oth->i; pj_oth = p_oth->j; 732968382fdSHong Zhang } else { 73320e3dc75SHong Zhang p_oth = NULL; 73420e3dc75SHong Zhang pi_oth = NULL; pj_oth = NULL; 735a7c7454dSHong Zhang } 736cf1a0672SHong Zhang 737cf1a0672SHong Zhang /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */ 738cf1a0672SHong Zhang /*-------------------------------------------------------------------*/ 739854ce69bSBarry Smith ierr = PetscMalloc1(am+2,&api);CHKERRQ(ierr); 740cf1a0672SHong Zhang ptap->api = api; 741cf1a0672SHong Zhang api[0] = 0; 742cf1a0672SHong Zhang 743cf1a0672SHong Zhang /* create and initialize a linked list */ 744cf1a0672SHong Zhang armax = ad->rmax+ao->rmax; 745a7c7454dSHong Zhang if (size >1) { 746cf1a0672SHong Zhang prmax = PetscMax(p_loc->rmax,p_oth->rmax); 747a7c7454dSHong Zhang } else { 748a7c7454dSHong Zhang prmax = p_loc->rmax; 749a7c7454dSHong Zhang } 750cf1a0672SHong Zhang nlnk_max = armax*prmax; 751cf1a0672SHong Zhang if (!nlnk_max || nlnk_max > pN) nlnk_max = pN; 752f84c536bSHong Zhang ierr = PetscLLCondensedCreate_Scalable(nlnk_max,&lnk);CHKERRQ(ierr); 753cf1a0672SHong Zhang 754cf1a0672SHong Zhang /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */ 755cf1a0672SHong Zhang ierr = PetscFreeSpaceGet((PetscInt)(fill*(adi[am]+aoi[am]+pi_loc[pm])),&free_space);CHKERRQ(ierr); 7562205254eSKarl Rupp 757cf1a0672SHong Zhang current_space = free_space; 758cf1a0672SHong Zhang 759cf1a0672SHong Zhang ierr = MatPreallocateInitialize(comm,am,pn,dnz,onz);CHKERRQ(ierr); 760cf1a0672SHong Zhang for (i=0; i<am; i++) { 761cf1a0672SHong Zhang /* diagonal portion of A */ 762cf1a0672SHong Zhang nzi = adi[i+1] - adi[i]; 763cf1a0672SHong Zhang for (j=0; j<nzi; j++) { 764cf1a0672SHong Zhang row = *adj++; 765cf1a0672SHong Zhang pnz = pi_loc[row+1] - pi_loc[row]; 766cf1a0672SHong Zhang Jptr = pj_loc + pi_loc[row]; 767cf1a0672SHong Zhang /* add non-zero cols of P into the sorted linked list lnk */ 768f84c536bSHong Zhang ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr); 769cf1a0672SHong Zhang } 770cf1a0672SHong Zhang /* off-diagonal portion of A */ 771cf1a0672SHong Zhang nzi = aoi[i+1] - aoi[i]; 772cf1a0672SHong Zhang for (j=0; j<nzi; j++) { 773cf1a0672SHong Zhang row = *aoj++; 774cf1a0672SHong Zhang pnz = pi_oth[row+1] - pi_oth[row]; 775cf1a0672SHong Zhang Jptr = pj_oth + pi_oth[row]; 776f84c536bSHong Zhang ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr); 777cf1a0672SHong Zhang } 778cf1a0672SHong Zhang 779f84c536bSHong Zhang apnz = *lnk; 780cf1a0672SHong Zhang api[i+1] = api[i] + apnz; 781cf1a0672SHong Zhang if (apnz > apnz_max) apnz_max = apnz; 782cf1a0672SHong Zhang 783cf1a0672SHong Zhang /* if free space is not available, double the total space in the list */ 784cf1a0672SHong Zhang if (current_space->local_remaining<apnz) { 785cf1a0672SHong Zhang ierr = PetscFreeSpaceGet(apnz+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 786cf1a0672SHong Zhang nspacedouble++; 787cf1a0672SHong Zhang } 788cf1a0672SHong Zhang 789cf1a0672SHong Zhang /* Copy data into free space, then initialize lnk */ 790f84c536bSHong Zhang ierr = PetscLLCondensedClean_Scalable(apnz,current_space->array,lnk);CHKERRQ(ierr); 791cf1a0672SHong Zhang ierr = MatPreallocateSet(i+rstart,apnz,current_space->array,dnz,onz);CHKERRQ(ierr); 7922205254eSKarl Rupp 793cf1a0672SHong Zhang current_space->array += apnz; 794cf1a0672SHong Zhang current_space->local_used += apnz; 795cf1a0672SHong Zhang current_space->local_remaining -= apnz; 796cf1a0672SHong Zhang } 797cf1a0672SHong Zhang 798cf1a0672SHong Zhang /* Allocate space for apj, initialize apj, and */ 799cf1a0672SHong Zhang /* destroy list of free space and other temporary array(s) */ 800854ce69bSBarry Smith ierr = PetscMalloc1(api[am]+1,&ptap->apj);CHKERRQ(ierr); 801cf1a0672SHong Zhang apj = ptap->apj; 802cf1a0672SHong Zhang ierr = PetscFreeSpaceContiguous(&free_space,ptap->apj);CHKERRQ(ierr); 803f84c536bSHong Zhang ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr); 804cf1a0672SHong Zhang 805cf1a0672SHong Zhang /* create and assemble symbolic parallel matrix Cmpi */ 806cf1a0672SHong Zhang /*----------------------------------------------------*/ 807cf1a0672SHong Zhang ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr); 808cf1a0672SHong Zhang ierr = MatSetSizes(Cmpi,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 80933d57670SJed Brown ierr = MatSetBlockSizesFromMats(Cmpi,A,P);CHKERRQ(ierr); 810cf1a0672SHong Zhang ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr); 811cf1a0672SHong Zhang ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr); 812cf1a0672SHong Zhang ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 813cf1a0672SHong Zhang 81429825010SHong Zhang /* malloc apa for assembly Cmpi */ 8151795a4d1SJed Brown ierr = PetscCalloc1(apnz_max,&apa);CHKERRQ(ierr); 8162205254eSKarl Rupp 81729825010SHong Zhang ptap->apa = apa; 818cf1a0672SHong Zhang for (i=0; i<am; i++) { 819cf1a0672SHong Zhang row = i + rstart; 820cf1a0672SHong Zhang apnz = api[i+1] - api[i]; 821cf1a0672SHong Zhang ierr = MatSetValues(Cmpi,1,&row,apnz,apj,apa,INSERT_VALUES);CHKERRQ(ierr); 822cf1a0672SHong Zhang apj += apnz; 823cf1a0672SHong Zhang } 824cf1a0672SHong Zhang ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 825cf1a0672SHong Zhang ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 826cf1a0672SHong Zhang 827cf1a0672SHong Zhang ptap->destroy = Cmpi->ops->destroy; 828cf1a0672SHong Zhang ptap->duplicate = Cmpi->ops->duplicate; 829f8efd71cSHong Zhang Cmpi->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ; 830cf1a0672SHong Zhang Cmpi->ops->destroy = MatDestroy_MPIAIJ_MatMatMult; 831cf1a0672SHong Zhang Cmpi->ops->duplicate = MatDuplicate_MPIAIJ_MatMatMult; 832cf1a0672SHong Zhang 833cf1a0672SHong Zhang /* attach the supporting struct to Cmpi for reuse */ 834cf1a0672SHong Zhang c = (Mat_MPIAIJ*)Cmpi->data; 835cf1a0672SHong Zhang c->ptap = ptap; 836cf1a0672SHong Zhang 837cf1a0672SHong Zhang *C = Cmpi; 838cf1a0672SHong Zhang 839cf1a0672SHong Zhang /* set MatInfo */ 840118727c9SMark F. Adams afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1) + 1.e-5; 841cf1a0672SHong Zhang if (afill < 1.0) afill = 1.0; 842cf1a0672SHong Zhang Cmpi->info.mallocs = nspacedouble; 843cf1a0672SHong Zhang Cmpi->info.fill_ratio_given = fill; 844cf1a0672SHong Zhang Cmpi->info.fill_ratio_needed = afill; 845cf1a0672SHong Zhang 846cf1a0672SHong Zhang #if defined(PETSC_USE_INFO) 847cf1a0672SHong Zhang if (api[am]) { 84857622a8eSBarry Smith ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr); 84957622a8eSBarry Smith ierr = PetscInfo1(Cmpi,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr); 850cf1a0672SHong Zhang } else { 851cf1a0672SHong Zhang ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr); 852cf1a0672SHong Zhang } 853cf1a0672SHong Zhang #endif 8541634b032SHong Zhang PetscFunctionReturn(0); 8551634b032SHong Zhang } 856f96d37fcSHong Zhang 857f96d37fcSHong Zhang /*-------------------------------------------------------------------------*/ 858f96d37fcSHong Zhang #undef __FUNCT__ 859f96d37fcSHong Zhang #define __FUNCT__ "MatTransposeMatMult_MPIAIJ_MPIAIJ" 860c5af039cSHong Zhang PetscErrorCode MatTransposeMatMult_MPIAIJ_MPIAIJ(Mat P,Mat A,MatReuse scall,PetscReal fill,Mat *C) 861f96d37fcSHong Zhang { 862f96d37fcSHong Zhang PetscErrorCode ierr; 863c216dbf3SHong Zhang const char *algTypes[3] = {"scalable","nonscalable","matmatmult"}; 864c216dbf3SHong Zhang PetscInt alg=0; /* set default algorithm */ 865f96d37fcSHong Zhang 866f96d37fcSHong Zhang PetscFunctionBegin; 867c216dbf3SHong Zhang if (scall == MAT_INITIAL_MATRIX) { 868c216dbf3SHong Zhang ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr); 869c216dbf3SHong Zhang ierr = PetscOptionsEList("-mattransposematmult_via","Algorithmic approach","MatTransposeMatMult",algTypes,3,algTypes[0],&alg,NULL);CHKERRQ(ierr); 870c216dbf3SHong Zhang ierr = PetscOptionsEnd();CHKERRQ(ierr); 871c216dbf3SHong Zhang 872c216dbf3SHong Zhang ierr = PetscLogEventBegin(MAT_TransposeMatMultSymbolic,P,A,0,0);CHKERRQ(ierr); 873c216dbf3SHong Zhang switch (alg) { 874c216dbf3SHong Zhang case 1: 8752bbb1c24SHong Zhang ierr = MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(P,A,fill,C);CHKERRQ(ierr); 876c216dbf3SHong Zhang break; 877c216dbf3SHong Zhang case 2: 878275476c6SMatthew G. Knepley { 879ab1b013aSHong Zhang Mat Pt; 880ab1b013aSHong Zhang Mat_PtAPMPI *ptap; 881ab1b013aSHong Zhang Mat_MPIAIJ *c; 882ab1b013aSHong Zhang ierr = MatTranspose(P,MAT_INITIAL_MATRIX,&Pt);CHKERRQ(ierr); 883ab1b013aSHong Zhang ierr = MatMatMult(Pt,A,MAT_INITIAL_MATRIX,fill,C);CHKERRQ(ierr); 884ab1b013aSHong Zhang c = (Mat_MPIAIJ*)(*C)->data; 885ab1b013aSHong Zhang ptap = c->ptap; 886ab1b013aSHong Zhang ptap->Pt = Pt; 887c216dbf3SHong Zhang (*C)->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult; 888c216dbf3SHong Zhang PetscFunctionReturn(0); 889275476c6SMatthew G. Knepley } 890c216dbf3SHong Zhang break; 891c216dbf3SHong Zhang default: 8926da69ca6SHong Zhang ierr = MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(P,A,fill,C);CHKERRQ(ierr); 893c216dbf3SHong Zhang break; 894c216dbf3SHong Zhang } 895c216dbf3SHong Zhang ierr = PetscLogEventEnd(MAT_TransposeMatMultSymbolic,P,A,0,0);CHKERRQ(ierr); 896c216dbf3SHong Zhang } 897c216dbf3SHong Zhang ierr = PetscLogEventBegin(MAT_TransposeMatMultNumeric,P,A,0,0);CHKERRQ(ierr); 898c216dbf3SHong Zhang ierr = (*(*C)->ops->mattransposemultnumeric)(P,A,*C);CHKERRQ(ierr); 899c216dbf3SHong Zhang ierr = PetscLogEventEnd(MAT_TransposeMatMultNumeric,P,A,0,0);CHKERRQ(ierr); 900ab1b013aSHong Zhang PetscFunctionReturn(0); 901ab1b013aSHong Zhang } 902ab1b013aSHong Zhang 903c216dbf3SHong Zhang /* This routine only works when scall=MAT_REUSE_MATRIX! */ 904c216dbf3SHong Zhang #undef __FUNCT__ 905c216dbf3SHong Zhang #define __FUNCT__ "MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult" 906c216dbf3SHong Zhang PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult(Mat P,Mat A,Mat C) 907c216dbf3SHong Zhang { 908c216dbf3SHong Zhang PetscErrorCode ierr; 9092bbb1c24SHong Zhang Mat_MPIAIJ *c=(Mat_MPIAIJ*)C->data; 9102bbb1c24SHong Zhang Mat_PtAPMPI *ptap= c->ptap; 9112bbb1c24SHong Zhang Mat Pt=ptap->Pt; 912c216dbf3SHong Zhang 913c216dbf3SHong Zhang PetscFunctionBegin; 914c216dbf3SHong Zhang ierr = MatTranspose(P,MAT_REUSE_MATRIX,&Pt);CHKERRQ(ierr); 915c216dbf3SHong Zhang ierr = MatMatMultNumeric(Pt,A,C);CHKERRQ(ierr); 916f96d37fcSHong Zhang PetscFunctionReturn(0); 917f96d37fcSHong Zhang } 918f96d37fcSHong Zhang 9196da69ca6SHong Zhang /* Non-scalable version, use dense axpy */ 920f96d37fcSHong Zhang #undef __FUNCT__ 9216da69ca6SHong Zhang #define __FUNCT__ "MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable" 9226da69ca6SHong Zhang PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable(Mat P,Mat A,Mat C) 923f96d37fcSHong Zhang { 924c5af039cSHong Zhang PetscErrorCode ierr; 925c5af039cSHong Zhang Mat_Merge_SeqsToMPI *merge; 926497f5370SHong Zhang Mat_MPIAIJ *p =(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data; 927c5af039cSHong Zhang Mat_SeqAIJ *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data; 928c5af039cSHong Zhang Mat_PtAPMPI *ptap; 929d6ab1dc0SHong Zhang PetscInt *adj,*aJ; 930497f5370SHong Zhang PetscInt i,j,k,anz,pnz,row,*cj; 931e745005bSHong Zhang MatScalar *ada,*aval,*ca,valtmp; 932c5af039cSHong Zhang PetscInt am =A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n; 933ce94432eSBarry Smith MPI_Comm comm; 934c5af039cSHong Zhang PetscMPIInt size,rank,taga,*len_s; 935c5af039cSHong Zhang PetscInt *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci; 936c5af039cSHong Zhang PetscInt **buf_ri,**buf_rj; 937c5af039cSHong Zhang PetscInt cnz=0,*bj_i,*bi,*bj,bnz,nextcj; /* bi,bj,ba: local array of C(mpi mat) */ 938c5af039cSHong Zhang MPI_Request *s_waits,*r_waits; 939c5af039cSHong Zhang MPI_Status *status; 940c5af039cSHong Zhang MatScalar **abuf_r,*ba_i,*pA,*coa,*ba; 941d6ab1dc0SHong Zhang PetscInt *ai,*aj,*coi,*coj; 942a2ea699eSBarry Smith PetscInt *poJ,*pdJ; 943c5af039cSHong Zhang Mat A_loc; 944c5af039cSHong Zhang Mat_SeqAIJ *a_loc; 945f96d37fcSHong Zhang 946f96d37fcSHong Zhang PetscFunctionBegin; 947ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr); 948c5af039cSHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 949c5af039cSHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 950c5af039cSHong Zhang 951c5af039cSHong Zhang ptap = c->ptap; 952c5af039cSHong Zhang merge = ptap->merge; 953c5af039cSHong Zhang 954c5af039cSHong Zhang /* 2) compute numeric C_seq = P_loc^T*A_loc*P - dominating part */ 955c5af039cSHong Zhang /*--------------------------------------------------------------*/ 956c5af039cSHong Zhang /* get data from symbolic products */ 957c5af039cSHong Zhang coi = merge->coi; coj = merge->coj; 958854ce69bSBarry Smith ierr = PetscCalloc1(coi[pon]+1,&coa);CHKERRQ(ierr); 959c5af039cSHong Zhang 960c5af039cSHong Zhang bi = merge->bi; bj = merge->bj; 961c5af039cSHong Zhang owners = merge->rowmap->range; 962854ce69bSBarry Smith ierr = PetscCalloc1(bi[cm]+1,&ba);CHKERRQ(ierr); 963c5af039cSHong Zhang 964c5af039cSHong Zhang /* get A_loc by taking all local rows of A */ 965c5af039cSHong Zhang A_loc = ptap->A_loc; 966c5af039cSHong Zhang ierr = MatMPIAIJGetLocalMat(A,MAT_REUSE_MATRIX,&A_loc);CHKERRQ(ierr); 967c5af039cSHong Zhang a_loc = (Mat_SeqAIJ*)(A_loc)->data; 968d6ab1dc0SHong Zhang ai = a_loc->i; 969d6ab1dc0SHong Zhang aj = a_loc->j; 970c5af039cSHong Zhang 971854ce69bSBarry Smith ierr = PetscCalloc1(A->cmap->N,&aval);CHKERRQ(ierr); /* non-scalable!!! */ 972c5af039cSHong Zhang 973c5af039cSHong Zhang for (i=0; i<am; i++) { 974e745005bSHong Zhang /* 2-a) put A[i,:] to dense array aval */ 975d6ab1dc0SHong Zhang anz = ai[i+1] - ai[i]; 976d6ab1dc0SHong Zhang adj = aj + ai[i]; 977d6ab1dc0SHong Zhang ada = a_loc->a + ai[i]; 978c5af039cSHong Zhang for (j=0; j<anz; j++) { 979e745005bSHong Zhang aval[adj[j]] = ada[j]; 980c5af039cSHong Zhang } 981c5af039cSHong Zhang 982c5af039cSHong Zhang /* 2-b) Compute Cseq = P_loc[i,:]^T*A[i,:] using outer product */ 983c5af039cSHong Zhang /*--------------------------------------------------------------*/ 984c5af039cSHong Zhang /* put the value into Co=(p->B)^T*A (off-diagonal part, send to others) */ 985c5af039cSHong Zhang pnz = po->i[i+1] - po->i[i]; 986c5af039cSHong Zhang poJ = po->j + po->i[i]; 987c5af039cSHong Zhang pA = po->a + po->i[i]; 988c5af039cSHong Zhang for (j=0; j<pnz; j++) { 989c5af039cSHong Zhang row = poJ[j]; 990c5af039cSHong Zhang cnz = coi[row+1] - coi[row]; 991c5af039cSHong Zhang cj = coj + coi[row]; 992c5af039cSHong Zhang ca = coa + coi[row]; 993c5af039cSHong Zhang /* perform dense axpy */ 994c5af039cSHong Zhang valtmp = pA[j]; 995c5af039cSHong Zhang for (k=0; k<cnz; k++) { 996e745005bSHong Zhang ca[k] += valtmp*aval[cj[k]]; 997c5af039cSHong Zhang } 998c5af039cSHong Zhang ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr); 999c5af039cSHong Zhang } 1000c5af039cSHong Zhang 1001c5af039cSHong Zhang /* put the value into Cd (diagonal part) */ 1002c5af039cSHong Zhang pnz = pd->i[i+1] - pd->i[i]; 1003c5af039cSHong Zhang pdJ = pd->j + pd->i[i]; 1004c5af039cSHong Zhang pA = pd->a + pd->i[i]; 1005c5af039cSHong Zhang for (j=0; j<pnz; j++) { 1006c5af039cSHong Zhang row = pdJ[j]; 1007c5af039cSHong Zhang cnz = bi[row+1] - bi[row]; 1008c5af039cSHong Zhang cj = bj + bi[row]; 1009c5af039cSHong Zhang ca = ba + bi[row]; 1010c5af039cSHong Zhang /* perform dense axpy */ 1011c5af039cSHong Zhang valtmp = pA[j]; 1012c5af039cSHong Zhang for (k=0; k<cnz; k++) { 1013e745005bSHong Zhang ca[k] += valtmp*aval[cj[k]]; 1014c5af039cSHong Zhang } 1015c5af039cSHong Zhang ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr); 1016c5af039cSHong Zhang } 1017c5af039cSHong Zhang 1018d6ab1dc0SHong Zhang /* zero the current row of Pt*A */ 1019d6ab1dc0SHong Zhang aJ = aj + ai[i]; 1020e745005bSHong Zhang for (k=0; k<anz; k++) aval[aJ[k]] = 0.0; 1021c5af039cSHong Zhang } 1022c5af039cSHong Zhang 1023c5af039cSHong Zhang /* 3) send and recv matrix values coa */ 1024c5af039cSHong Zhang /*------------------------------------*/ 1025c5af039cSHong Zhang buf_ri = merge->buf_ri; 1026c5af039cSHong Zhang buf_rj = merge->buf_rj; 1027c5af039cSHong Zhang len_s = merge->len_s; 1028c5af039cSHong Zhang ierr = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr); 1029c5af039cSHong Zhang ierr = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr); 1030c5af039cSHong Zhang 1031dcca6d9dSJed Brown ierr = PetscMalloc2(merge->nsend+1,&s_waits,size,&status);CHKERRQ(ierr); 1032c5af039cSHong Zhang for (proc=0,k=0; proc<size; proc++) { 1033c5af039cSHong Zhang if (!len_s[proc]) continue; 1034c5af039cSHong Zhang i = merge->owners_co[proc]; 1035c5af039cSHong Zhang ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr); 1036c5af039cSHong Zhang k++; 1037c5af039cSHong Zhang } 1038c5af039cSHong Zhang if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);} 1039c5af039cSHong Zhang if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);} 1040c5af039cSHong Zhang 1041c5af039cSHong Zhang ierr = PetscFree2(s_waits,status);CHKERRQ(ierr); 1042c5af039cSHong Zhang ierr = PetscFree(r_waits);CHKERRQ(ierr); 1043c5af039cSHong Zhang ierr = PetscFree(coa);CHKERRQ(ierr); 1044c5af039cSHong Zhang 1045c5af039cSHong Zhang /* 4) insert local Cseq and received values into Cmpi */ 1046c5af039cSHong Zhang /*----------------------------------------------------*/ 1047dcca6d9dSJed Brown ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr); 1048c5af039cSHong Zhang for (k=0; k<merge->nrecv; k++) { 1049c36aecf5SHong Zhang buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 1050c5af039cSHong Zhang nrows = *(buf_ri_k[k]); 1051c5af039cSHong Zhang nextrow[k] = buf_ri_k[k]+1; /* next row number of k-th recved i-structure */ 1052c5af039cSHong Zhang nextci[k] = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure */ 1053c5af039cSHong Zhang } 1054c5af039cSHong Zhang 1055c5af039cSHong Zhang for (i=0; i<cm; i++) { 1056c5af039cSHong Zhang row = owners[rank] + i; /* global row index of C_seq */ 1057c5af039cSHong Zhang bj_i = bj + bi[i]; /* col indices of the i-th row of C */ 1058c5af039cSHong Zhang ba_i = ba + bi[i]; 1059c5af039cSHong Zhang bnz = bi[i+1] - bi[i]; 1060c5af039cSHong Zhang /* add received vals into ba */ 1061c5af039cSHong Zhang for (k=0; k<merge->nrecv; k++) { /* k-th received message */ 1062c5af039cSHong Zhang /* i-th row */ 1063c5af039cSHong Zhang if (i == *nextrow[k]) { 1064c5af039cSHong Zhang cnz = *(nextci[k]+1) - *nextci[k]; 1065c5af039cSHong Zhang cj = buf_rj[k] + *(nextci[k]); 1066c5af039cSHong Zhang ca = abuf_r[k] + *(nextci[k]); 1067c5af039cSHong Zhang nextcj = 0; 1068c5af039cSHong Zhang for (j=0; nextcj<cnz; j++) { 1069c5af039cSHong Zhang if (bj_i[j] == cj[nextcj]) { /* bcol == ccol */ 1070c5af039cSHong Zhang ba_i[j] += ca[nextcj++]; 1071c5af039cSHong Zhang } 1072c5af039cSHong Zhang } 1073c5af039cSHong Zhang nextrow[k]++; nextci[k]++; 1074c5af039cSHong Zhang ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr); 1075c5af039cSHong Zhang } 1076c5af039cSHong Zhang } 1077c5af039cSHong Zhang ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr); 1078c5af039cSHong Zhang } 1079c5af039cSHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1080c5af039cSHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1081c5af039cSHong Zhang 1082c5af039cSHong Zhang ierr = PetscFree(ba);CHKERRQ(ierr); 1083c5af039cSHong Zhang ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr); 1084c5af039cSHong Zhang ierr = PetscFree(abuf_r);CHKERRQ(ierr); 1085c5af039cSHong Zhang ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr); 1086e745005bSHong Zhang ierr = PetscFree(aval);CHKERRQ(ierr); 1087f96d37fcSHong Zhang PetscFunctionReturn(0); 1088f96d37fcSHong Zhang } 1089f96d37fcSHong Zhang 1090c9ba551fSBarry Smith PetscErrorCode MatDuplicate_MPIAIJ_MatPtAP(Mat, MatDuplicateOption,Mat*); 1091c5af039cSHong Zhang /* This routine is modified from MatPtAPSymbolic_MPIAIJ_MPIAIJ() */ 1092f96d37fcSHong Zhang #undef __FUNCT__ 10932bbb1c24SHong Zhang #define __FUNCT__ "MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable" 10942bbb1c24SHong Zhang PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat P,Mat A,PetscReal fill,Mat *C) 1095f96d37fcSHong Zhang { 1096f96d37fcSHong Zhang PetscErrorCode ierr; 10974e938277SHong Zhang Mat Cmpi,A_loc,POt,PDt; 1098f96d37fcSHong Zhang Mat_PtAPMPI *ptap; 10990298fd71SBarry Smith PetscFreeSpaceList free_space=NULL,current_space=NULL; 1100f96d37fcSHong Zhang Mat_MPIAIJ *p =(Mat_MPIAIJ*)P->data,*c; 1101f96d37fcSHong Zhang PetscInt *pdti,*pdtj,*poti,*potj,*ptJ; 1102f96d37fcSHong Zhang PetscInt nnz; 11034e938277SHong Zhang PetscInt *lnk,*owners_co,*coi,*coj,i,k,pnz,row; 1104497f5370SHong Zhang PetscInt am=A->rmap->n,pn=P->cmap->n; 1105f96d37fcSHong Zhang PetscBT lnkbt; 1106ce94432eSBarry Smith MPI_Comm comm; 1107f96d37fcSHong Zhang PetscMPIInt size,rank,tagi,tagj,*len_si,*len_s,*len_ri; 1108f96d37fcSHong Zhang PetscInt **buf_rj,**buf_ri,**buf_ri_k; 1109f96d37fcSHong Zhang PetscInt len,proc,*dnz,*onz,*owners; 1110f96d37fcSHong Zhang PetscInt nzi,*bi,*bj; 1111f96d37fcSHong Zhang PetscInt nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci; 1112f96d37fcSHong Zhang MPI_Request *swaits,*rwaits; 1113f96d37fcSHong Zhang MPI_Status *sstatus,rstatus; 1114f96d37fcSHong Zhang Mat_Merge_SeqsToMPI *merge; 1115d6ab1dc0SHong Zhang PetscInt *ai,*aj,*Jptr,anz,*prmap=p->garray,pon,nspacedouble=0,j; 1116f96d37fcSHong Zhang PetscReal afill =1.0,afill_tmp; 1117438d860cSHong Zhang PetscInt rstart = P->cmap->rstart,rmax,aN=A->cmap->N; 1118f96d37fcSHong Zhang PetscScalar *vals; 11194e938277SHong Zhang Mat_SeqAIJ *a_loc, *pdt,*pot; 1120f96d37fcSHong Zhang 1121f96d37fcSHong Zhang PetscFunctionBegin; 1122ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 1123c5af039cSHong Zhang /* check if matrix local sizes are compatible */ 1124c5af039cSHong Zhang if (A->rmap->rstart != P->rmap->rstart || A->rmap->rend != P->rmap->rend) { 1125c5af039cSHong 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); 1126c5af039cSHong Zhang } 1127c5af039cSHong Zhang 1128f96d37fcSHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1129f96d37fcSHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 1130f96d37fcSHong Zhang 1131f96d37fcSHong Zhang /* create struct Mat_PtAPMPI and attached it to C later */ 1132b00a9115SJed Brown ierr = PetscNew(&ptap);CHKERRQ(ierr); 1133f96d37fcSHong Zhang 1134f96d37fcSHong Zhang /* get A_loc by taking all local rows of A */ 1135f96d37fcSHong Zhang ierr = MatMPIAIJGetLocalMat(A,MAT_INITIAL_MATRIX,&A_loc);CHKERRQ(ierr); 11362205254eSKarl Rupp 1137c5af039cSHong Zhang ptap->A_loc = A_loc; 11382205254eSKarl Rupp 11391c7d5954SHong Zhang a_loc = (Mat_SeqAIJ*)(A_loc)->data; 1140d6ab1dc0SHong Zhang ai = a_loc->i; 1141d6ab1dc0SHong Zhang aj = a_loc->j; 1142f96d37fcSHong Zhang 1143f96d37fcSHong Zhang /* determine symbolic Co=(p->B)^T*A - send to others */ 1144f96d37fcSHong Zhang /*----------------------------------------------------*/ 11454e938277SHong Zhang ierr = MatTransposeSymbolic_SeqAIJ(p->A,&PDt);CHKERRQ(ierr); 11464e938277SHong Zhang pdt = (Mat_SeqAIJ*)PDt->data; 11474e938277SHong Zhang pdti = pdt->i; pdtj = pdt->j; 11484e938277SHong Zhang 11494e938277SHong Zhang ierr = MatTransposeSymbolic_SeqAIJ(p->B,&POt);CHKERRQ(ierr); 11504e938277SHong Zhang pot = (Mat_SeqAIJ*)POt->data; 11514e938277SHong Zhang poti = pot->i; potj = pot->j; 1152f96d37fcSHong Zhang 1153f96d37fcSHong Zhang /* then, compute symbolic Co = (p->B)^T*A */ 11542205254eSKarl Rupp pon = (p->B)->cmap->n; /* total num of rows to be sent to other processors >= (num of nonzero rows of C_seq) - pn */ 1155854ce69bSBarry Smith ierr = PetscMalloc1(pon+1,&coi);CHKERRQ(ierr); 1156f96d37fcSHong Zhang coi[0] = 0; 1157f96d37fcSHong Zhang 1158f96d37fcSHong Zhang /* set initial free space to be fill*(nnz(p->B) + nnz(A)) */ 1159d6ab1dc0SHong Zhang nnz = fill*(poti[pon] + ai[am]); 1160a2ea699eSBarry Smith ierr = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr); 1161f96d37fcSHong Zhang current_space = free_space; 116219f0e57aSHong Zhang 1163f96d37fcSHong Zhang /* create and initialize a linked list */ 1164438d860cSHong Zhang ierr = PetscLLCondensedCreate(aN,aN,&lnk,&lnkbt);CHKERRQ(ierr); 1165f96d37fcSHong Zhang 1166f96d37fcSHong Zhang for (i=0; i<pon; i++) { 1167f96d37fcSHong Zhang pnz = poti[i+1] - poti[i]; 1168f96d37fcSHong Zhang ptJ = potj + poti[i]; 1169f96d37fcSHong Zhang for (j=0; j<pnz; j++) { 1170f96d37fcSHong Zhang row = ptJ[j]; /* row of A_loc == col of Pot */ 1171d6ab1dc0SHong Zhang anz = ai[row+1] - ai[row]; 1172d6ab1dc0SHong Zhang Jptr = aj + ai[row]; 1173f96d37fcSHong Zhang /* add non-zero cols of AP into the sorted linked list lnk */ 1174d6ab1dc0SHong Zhang ierr = PetscLLCondensedAddSorted(anz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 1175f96d37fcSHong Zhang } 11764e938277SHong Zhang nnz = lnk[0]; 1177f96d37fcSHong Zhang 1178f96d37fcSHong Zhang /* If free space is not available, double the total space in the list */ 1179f96d37fcSHong Zhang if (current_space->local_remaining<nnz) { 1180f96d37fcSHong Zhang ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 1181f96d37fcSHong Zhang nspacedouble++; 1182f96d37fcSHong Zhang } 1183f96d37fcSHong Zhang 1184f96d37fcSHong Zhang /* Copy data into free space, and zero out denserows */ 11854e938277SHong Zhang ierr = PetscLLCondensedClean(aN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 11862205254eSKarl Rupp 1187f96d37fcSHong Zhang current_space->array += nnz; 1188f96d37fcSHong Zhang current_space->local_used += nnz; 1189f96d37fcSHong Zhang current_space->local_remaining -= nnz; 11902205254eSKarl Rupp 1191f96d37fcSHong Zhang coi[i+1] = coi[i] + nnz; 1192f96d37fcSHong Zhang } 1193f96d37fcSHong Zhang 1194854ce69bSBarry Smith ierr = PetscMalloc1(coi[pon]+1,&coj);CHKERRQ(ierr); 1195f96d37fcSHong Zhang ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr); 11962205254eSKarl Rupp 1197118727c9SMark F. Adams afill_tmp = (PetscReal)coi[pon]/(poti[pon] + ai[am]+1); 1198f96d37fcSHong Zhang if (afill_tmp > afill) afill = afill_tmp; 1199f96d37fcSHong Zhang 1200f96d37fcSHong Zhang /* send j-array (coj) of Co to other processors */ 1201f96d37fcSHong Zhang /*----------------------------------------------*/ 1202f96d37fcSHong Zhang /* determine row ownership */ 1203b00a9115SJed Brown ierr = PetscNew(&merge);CHKERRQ(ierr); 1204f96d37fcSHong Zhang ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr); 12052205254eSKarl Rupp 1206f96d37fcSHong Zhang merge->rowmap->n = pn; 1207f96d37fcSHong Zhang merge->rowmap->bs = 1; 12082205254eSKarl Rupp 1209f96d37fcSHong Zhang ierr = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr); 1210f96d37fcSHong Zhang owners = merge->rowmap->range; 1211f96d37fcSHong Zhang 1212f96d37fcSHong Zhang /* determine the number of messages to send, their lengths */ 12131795a4d1SJed Brown ierr = PetscCalloc1(size,&len_si);CHKERRQ(ierr); 1214785e854fSJed Brown ierr = PetscMalloc1(size,&merge->len_s);CHKERRQ(ierr); 12152205254eSKarl Rupp 1216f96d37fcSHong Zhang len_s = merge->len_s; 1217f96d37fcSHong Zhang merge->nsend = 0; 1218f96d37fcSHong Zhang 1219854ce69bSBarry Smith ierr = PetscMalloc1(size+2,&owners_co);CHKERRQ(ierr); 1220f96d37fcSHong Zhang ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 1221f96d37fcSHong Zhang 1222f96d37fcSHong Zhang proc = 0; 1223f96d37fcSHong Zhang for (i=0; i<pon; i++) { 1224f96d37fcSHong Zhang while (prmap[i] >= owners[proc+1]) proc++; 1225f96d37fcSHong Zhang len_si[proc]++; /* num of rows in Co to be sent to [proc] */ 1226f96d37fcSHong Zhang len_s[proc] += coi[i+1] - coi[i]; 1227f96d37fcSHong Zhang } 1228f96d37fcSHong Zhang 1229f96d37fcSHong Zhang len = 0; /* max length of buf_si[] */ 1230f96d37fcSHong Zhang owners_co[0] = 0; 1231f96d37fcSHong Zhang for (proc=0; proc<size; proc++) { 1232f96d37fcSHong Zhang owners_co[proc+1] = owners_co[proc] + len_si[proc]; 1233f96d37fcSHong Zhang if (len_si[proc]) { 1234f96d37fcSHong Zhang merge->nsend++; 1235f96d37fcSHong Zhang len_si[proc] = 2*(len_si[proc] + 1); 1236f96d37fcSHong Zhang len += len_si[proc]; 1237f96d37fcSHong Zhang } 1238f96d37fcSHong Zhang } 1239f96d37fcSHong Zhang 1240f96d37fcSHong Zhang /* determine the number and length of messages to receive for coi and coj */ 12410298fd71SBarry Smith ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&merge->nrecv);CHKERRQ(ierr); 1242f96d37fcSHong Zhang ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr); 1243f96d37fcSHong Zhang 1244f96d37fcSHong Zhang /* post the Irecv and Isend of coj */ 1245f96d37fcSHong Zhang ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr); 1246f96d37fcSHong Zhang ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr); 1247854ce69bSBarry Smith ierr = PetscMalloc1(merge->nsend+1,&swaits);CHKERRQ(ierr); 1248f96d37fcSHong Zhang for (proc=0, k=0; proc<size; proc++) { 1249f96d37fcSHong Zhang if (!len_s[proc]) continue; 1250f96d37fcSHong Zhang i = owners_co[proc]; 1251f96d37fcSHong Zhang ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr); 1252f96d37fcSHong Zhang k++; 1253f96d37fcSHong Zhang } 1254f96d37fcSHong Zhang 1255f96d37fcSHong Zhang /* receives and sends of coj are complete */ 1256785e854fSJed Brown ierr = PetscMalloc1(size,&sstatus);CHKERRQ(ierr); 1257f96d37fcSHong Zhang for (i=0; i<merge->nrecv; i++) { 1258c280ed6aSJed Brown PetscMPIInt icompleted; 1259c280ed6aSJed Brown ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 1260f96d37fcSHong Zhang } 1261f96d37fcSHong Zhang ierr = PetscFree(rwaits);CHKERRQ(ierr); 1262f96d37fcSHong Zhang if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);} 1263f96d37fcSHong Zhang 1264f96d37fcSHong Zhang /* send and recv coi */ 1265f96d37fcSHong Zhang /*-------------------*/ 1266f96d37fcSHong Zhang ierr = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr); 1267f96d37fcSHong Zhang ierr = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr); 1268854ce69bSBarry Smith ierr = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr); 1269f96d37fcSHong Zhang buf_si = buf_s; /* points to the beginning of k-th msg to be sent */ 1270f96d37fcSHong Zhang for (proc=0,k=0; proc<size; proc++) { 1271f96d37fcSHong Zhang if (!len_s[proc]) continue; 1272f96d37fcSHong Zhang /* form outgoing message for i-structure: 1273f96d37fcSHong Zhang buf_si[0]: nrows to be sent 1274f96d37fcSHong Zhang [1:nrows]: row index (global) 1275f96d37fcSHong Zhang [nrows+1:2*nrows+1]: i-structure index 1276f96d37fcSHong Zhang */ 1277f96d37fcSHong Zhang /*-------------------------------------------*/ 1278f96d37fcSHong Zhang nrows = len_si[proc]/2 - 1; 1279f96d37fcSHong Zhang buf_si_i = buf_si + nrows+1; 1280f96d37fcSHong Zhang buf_si[0] = nrows; 1281f96d37fcSHong Zhang buf_si_i[0] = 0; 1282f96d37fcSHong Zhang nrows = 0; 1283f96d37fcSHong Zhang for (i=owners_co[proc]; i<owners_co[proc+1]; i++) { 1284f96d37fcSHong Zhang nzi = coi[i+1] - coi[i]; 1285f96d37fcSHong Zhang buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */ 1286f96d37fcSHong Zhang buf_si[nrows+1] = prmap[i] -owners[proc]; /* local row index */ 1287f96d37fcSHong Zhang nrows++; 1288f96d37fcSHong Zhang } 1289f96d37fcSHong Zhang ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr); 1290f96d37fcSHong Zhang k++; 1291f96d37fcSHong Zhang buf_si += len_si[proc]; 1292f96d37fcSHong Zhang } 1293f96d37fcSHong Zhang i = merge->nrecv; 1294f96d37fcSHong Zhang while (i--) { 1295c280ed6aSJed Brown PetscMPIInt icompleted; 1296c280ed6aSJed Brown ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 1297f96d37fcSHong Zhang } 1298f96d37fcSHong Zhang ierr = PetscFree(rwaits);CHKERRQ(ierr); 1299f96d37fcSHong Zhang if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);} 1300f96d37fcSHong Zhang ierr = PetscFree(len_si);CHKERRQ(ierr); 1301f96d37fcSHong Zhang ierr = PetscFree(len_ri);CHKERRQ(ierr); 1302f96d37fcSHong Zhang ierr = PetscFree(swaits);CHKERRQ(ierr); 1303f96d37fcSHong Zhang ierr = PetscFree(sstatus);CHKERRQ(ierr); 1304f96d37fcSHong Zhang ierr = PetscFree(buf_s);CHKERRQ(ierr); 1305f96d37fcSHong Zhang 1306f96d37fcSHong Zhang /* compute the local portion of C (mpi mat) */ 1307f96d37fcSHong Zhang /*------------------------------------------*/ 1308f96d37fcSHong Zhang /* allocate bi array and free space for accumulating nonzero column info */ 1309854ce69bSBarry Smith ierr = PetscMalloc1(pn+1,&bi);CHKERRQ(ierr); 1310f96d37fcSHong Zhang bi[0] = 0; 1311f96d37fcSHong Zhang 1312c36aecf5SHong Zhang /* set initial free space to be fill*(nnz(P) + nnz(A)) */ 1313d6ab1dc0SHong Zhang nnz = fill*(pdti[pn] + poti[pon] + ai[am]); 1314a2ea699eSBarry Smith ierr = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr); 1315f96d37fcSHong Zhang current_space = free_space; 1316f96d37fcSHong Zhang 1317dcca6d9dSJed Brown ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr); 1318f96d37fcSHong Zhang for (k=0; k<merge->nrecv; k++) { 1319f96d37fcSHong Zhang buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 1320f96d37fcSHong Zhang nrows = *buf_ri_k[k]; 1321f96d37fcSHong Zhang nextrow[k] = buf_ri_k[k] + 1; /* next row number of k-th recved i-structure */ 1322f96d37fcSHong Zhang nextci[k] = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure */ 1323f96d37fcSHong Zhang } 1324f96d37fcSHong Zhang 13251c7d5954SHong Zhang ierr = MatPreallocateInitialize(comm,pn,A->cmap->n,dnz,onz);CHKERRQ(ierr); 1326f96d37fcSHong Zhang rmax = 0; 1327f96d37fcSHong Zhang for (i=0; i<pn; i++) { 1328f96d37fcSHong Zhang /* add pdt[i,:]*AP into lnk */ 1329f96d37fcSHong Zhang pnz = pdti[i+1] - pdti[i]; 1330f96d37fcSHong Zhang ptJ = pdtj + pdti[i]; 1331f96d37fcSHong Zhang for (j=0; j<pnz; j++) { 1332f96d37fcSHong Zhang row = ptJ[j]; /* row of AP == col of Pt */ 1333d6ab1dc0SHong Zhang anz = ai[row+1] - ai[row]; 1334d6ab1dc0SHong Zhang Jptr = aj + ai[row]; 1335f96d37fcSHong Zhang /* add non-zero cols of AP into the sorted linked list lnk */ 1336d6ab1dc0SHong Zhang ierr = PetscLLCondensedAddSorted(anz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 1337f96d37fcSHong Zhang } 13384e938277SHong Zhang 1339f96d37fcSHong Zhang /* add received col data into lnk */ 1340f96d37fcSHong Zhang for (k=0; k<merge->nrecv; k++) { /* k-th received message */ 1341f96d37fcSHong Zhang if (i == *nextrow[k]) { /* i-th row */ 1342f96d37fcSHong Zhang nzi = *(nextci[k]+1) - *nextci[k]; 1343f96d37fcSHong Zhang Jptr = buf_rj[k] + *nextci[k]; 13444e938277SHong Zhang ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr); 1345f96d37fcSHong Zhang nextrow[k]++; nextci[k]++; 1346f96d37fcSHong Zhang } 1347f96d37fcSHong Zhang } 13484e938277SHong Zhang nnz = lnk[0]; 1349f96d37fcSHong Zhang 1350f96d37fcSHong Zhang /* if free space is not available, make more free space */ 1351f96d37fcSHong Zhang if (current_space->local_remaining<nnz) { 1352f96d37fcSHong Zhang ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 1353f96d37fcSHong Zhang nspacedouble++; 1354f96d37fcSHong Zhang } 1355f96d37fcSHong Zhang /* copy data into free space, then initialize lnk */ 13564e938277SHong Zhang ierr = PetscLLCondensedClean(aN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 1357f96d37fcSHong Zhang ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr); 13582205254eSKarl Rupp 1359f96d37fcSHong Zhang current_space->array += nnz; 1360f96d37fcSHong Zhang current_space->local_used += nnz; 1361f96d37fcSHong Zhang current_space->local_remaining -= nnz; 13622205254eSKarl Rupp 1363f96d37fcSHong Zhang bi[i+1] = bi[i] + nnz; 1364f96d37fcSHong Zhang if (nnz > rmax) rmax = nnz; 1365f96d37fcSHong Zhang } 1366f96d37fcSHong Zhang ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr); 1367f96d37fcSHong Zhang 1368854ce69bSBarry Smith ierr = PetscMalloc1(bi[pn]+1,&bj);CHKERRQ(ierr); 1369f96d37fcSHong Zhang ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr); 13702205254eSKarl Rupp 1371118727c9SMark F. Adams afill_tmp = (PetscReal)bi[pn]/(pdti[pn] + poti[pon] + ai[am]+1); 1372f96d37fcSHong Zhang if (afill_tmp > afill) afill = afill_tmp; 1373d6ab1dc0SHong Zhang ierr = PetscLLCondensedDestroy(lnk,lnkbt);CHKERRQ(ierr); 13744e938277SHong Zhang ierr = MatDestroy(&POt);CHKERRQ(ierr); 13754e938277SHong Zhang ierr = MatDestroy(&PDt);CHKERRQ(ierr); 1376f96d37fcSHong Zhang 13771c7d5954SHong Zhang /* create symbolic parallel matrix Cmpi - why cannot be assembled in Numeric part */ 13781c7d5954SHong Zhang /*----------------------------------------------------------------------------------*/ 13791795a4d1SJed Brown ierr = PetscCalloc1(rmax+1,&vals);CHKERRQ(ierr); 1380f96d37fcSHong Zhang 1381f96d37fcSHong Zhang ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr); 1382f96d37fcSHong Zhang ierr = MatSetSizes(Cmpi,pn,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 138333d57670SJed Brown ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(A->cmap->bs));CHKERRQ(ierr); 1384f96d37fcSHong Zhang ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr); 1385f96d37fcSHong Zhang ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr); 1386f96d37fcSHong Zhang ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 1387f96d37fcSHong Zhang ierr = MatSetBlockSize(Cmpi,1);CHKERRQ(ierr); 1388f96d37fcSHong Zhang for (i=0; i<pn; i++) { 1389f96d37fcSHong Zhang row = i + rstart; 1390f96d37fcSHong Zhang nnz = bi[i+1] - bi[i]; 1391f96d37fcSHong Zhang Jptr = bj + bi[i]; 1392f96d37fcSHong Zhang ierr = MatSetValues(Cmpi,1,&row,nnz,Jptr,vals,INSERT_VALUES);CHKERRQ(ierr); 1393f96d37fcSHong Zhang } 1394f96d37fcSHong Zhang ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1395f96d37fcSHong Zhang ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1396f96d37fcSHong Zhang ierr = PetscFree(vals);CHKERRQ(ierr); 1397f96d37fcSHong Zhang 1398f96d37fcSHong Zhang merge->bi = bi; 1399f96d37fcSHong Zhang merge->bj = bj; 1400f96d37fcSHong Zhang merge->coi = coi; 1401f96d37fcSHong Zhang merge->coj = coj; 1402f96d37fcSHong Zhang merge->buf_ri = buf_ri; 1403f96d37fcSHong Zhang merge->buf_rj = buf_rj; 1404f96d37fcSHong Zhang merge->owners_co = owners_co; 1405f96d37fcSHong Zhang merge->destroy = Cmpi->ops->destroy; 1406f96d37fcSHong Zhang merge->duplicate = Cmpi->ops->duplicate; 1407f96d37fcSHong Zhang 14086da69ca6SHong Zhang Cmpi->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable; 1409c5af039cSHong Zhang Cmpi->ops->destroy = MatDestroy_MPIAIJ_PtAP; 1410c9ba551fSBarry Smith Cmpi->ops->duplicate = MatDuplicate_MPIAIJ_MatPtAP; 1411f96d37fcSHong Zhang 1412f96d37fcSHong Zhang /* attach the supporting struct to Cmpi for reuse */ 1413f96d37fcSHong Zhang c = (Mat_MPIAIJ*)Cmpi->data; 1414f96d37fcSHong Zhang c->ptap = ptap; 14150298fd71SBarry Smith ptap->api = NULL; 14160298fd71SBarry Smith ptap->apj = NULL; 1417f96d37fcSHong Zhang ptap->merge = merge; 1418d6ab1dc0SHong Zhang ptap->rmax = rmax; 1419d6ab1dc0SHong Zhang 1420d6ab1dc0SHong Zhang *C = Cmpi; 1421d6ab1dc0SHong Zhang #if defined(PETSC_USE_INFO) 1422d6ab1dc0SHong Zhang if (bi[pn] != 0) { 142357622a8eSBarry Smith ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr); 142457622a8eSBarry Smith ierr = PetscInfo1(Cmpi,"Use MatTransposeMatMult(A,B,MatReuse,%g,&C) for best performance.\n",(double)afill);CHKERRQ(ierr); 1425d6ab1dc0SHong Zhang } else { 1426d6ab1dc0SHong Zhang ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr); 1427d6ab1dc0SHong Zhang } 1428d6ab1dc0SHong Zhang #endif 1429d6ab1dc0SHong Zhang PetscFunctionReturn(0); 1430d6ab1dc0SHong Zhang } 1431d6ab1dc0SHong Zhang 1432d6ab1dc0SHong Zhang #undef __FUNCT__ 14336da69ca6SHong Zhang #define __FUNCT__ "MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ" 14346da69ca6SHong Zhang PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ(Mat P,Mat A,Mat C) 1435d6ab1dc0SHong Zhang { 1436d6ab1dc0SHong Zhang PetscErrorCode ierr; 1437d6ab1dc0SHong Zhang Mat_Merge_SeqsToMPI *merge; 1438d6ab1dc0SHong Zhang Mat_MPIAIJ *p =(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data; 1439d6ab1dc0SHong Zhang Mat_SeqAIJ *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data; 1440d6ab1dc0SHong Zhang Mat_PtAPMPI *ptap; 1441e745005bSHong Zhang PetscInt *adj; 1442e745005bSHong Zhang PetscInt i,j,k,anz,pnz,row,*cj,nexta; 1443e745005bSHong Zhang MatScalar *ada,*ca,valtmp; 1444d6ab1dc0SHong Zhang PetscInt am =A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n; 1445ce94432eSBarry Smith MPI_Comm comm; 1446d6ab1dc0SHong Zhang PetscMPIInt size,rank,taga,*len_s; 1447d6ab1dc0SHong Zhang PetscInt *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci; 1448d6ab1dc0SHong Zhang PetscInt **buf_ri,**buf_rj; 1449d6ab1dc0SHong Zhang PetscInt cnz=0,*bj_i,*bi,*bj,bnz,nextcj; /* bi,bj,ba: local array of C(mpi mat) */ 1450d6ab1dc0SHong Zhang MPI_Request *s_waits,*r_waits; 1451d6ab1dc0SHong Zhang MPI_Status *status; 1452d6ab1dc0SHong Zhang MatScalar **abuf_r,*ba_i,*pA,*coa,*ba; 1453d6ab1dc0SHong Zhang PetscInt *ai,*aj,*coi,*coj; 1454a2ea699eSBarry Smith PetscInt *poJ,*pdJ; 1455d6ab1dc0SHong Zhang Mat A_loc; 1456d6ab1dc0SHong Zhang Mat_SeqAIJ *a_loc; 1457d6ab1dc0SHong Zhang 1458d6ab1dc0SHong Zhang PetscFunctionBegin; 1459ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr); 1460d6ab1dc0SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1461d6ab1dc0SHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 1462d6ab1dc0SHong Zhang 1463d6ab1dc0SHong Zhang ptap = c->ptap; 1464d6ab1dc0SHong Zhang merge = ptap->merge; 1465d6ab1dc0SHong Zhang 1466e745005bSHong Zhang /* 2) compute numeric C_seq = P_loc^T*A_loc */ 1467e745005bSHong Zhang /*------------------------------------------*/ 1468d6ab1dc0SHong Zhang /* get data from symbolic products */ 1469d6ab1dc0SHong Zhang coi = merge->coi; coj = merge->coj; 1470854ce69bSBarry Smith ierr = PetscCalloc1(coi[pon]+1,&coa);CHKERRQ(ierr); 1471d6ab1dc0SHong Zhang bi = merge->bi; bj = merge->bj; 1472d6ab1dc0SHong Zhang owners = merge->rowmap->range; 14731795a4d1SJed Brown ierr = PetscCalloc1(bi[cm]+1,&ba);CHKERRQ(ierr); 1474d6ab1dc0SHong Zhang 1475d6ab1dc0SHong Zhang /* get A_loc by taking all local rows of A */ 1476d6ab1dc0SHong Zhang A_loc = ptap->A_loc; 1477d6ab1dc0SHong Zhang ierr = MatMPIAIJGetLocalMat(A,MAT_REUSE_MATRIX,&A_loc);CHKERRQ(ierr); 1478d6ab1dc0SHong Zhang a_loc = (Mat_SeqAIJ*)(A_loc)->data; 1479d6ab1dc0SHong Zhang ai = a_loc->i; 1480d6ab1dc0SHong Zhang aj = a_loc->j; 1481d6ab1dc0SHong Zhang 1482d6ab1dc0SHong Zhang for (i=0; i<am; i++) { 1483d6ab1dc0SHong Zhang anz = ai[i+1] - ai[i]; 1484d6ab1dc0SHong Zhang adj = aj + ai[i]; 1485d6ab1dc0SHong Zhang ada = a_loc->a + ai[i]; 1486d6ab1dc0SHong Zhang 1487d6ab1dc0SHong Zhang /* 2-b) Compute Cseq = P_loc[i,:]^T*A[i,:] using outer product */ 1488e745005bSHong Zhang /*-------------------------------------------------------------*/ 1489d6ab1dc0SHong Zhang /* put the value into Co=(p->B)^T*A (off-diagonal part, send to others) */ 1490d6ab1dc0SHong Zhang pnz = po->i[i+1] - po->i[i]; 1491d6ab1dc0SHong Zhang poJ = po->j + po->i[i]; 1492d6ab1dc0SHong Zhang pA = po->a + po->i[i]; 1493d6ab1dc0SHong Zhang for (j=0; j<pnz; j++) { 1494d6ab1dc0SHong Zhang row = poJ[j]; 1495d6ab1dc0SHong Zhang cj = coj + coi[row]; 1496d6ab1dc0SHong Zhang ca = coa + coi[row]; 1497e745005bSHong Zhang /* perform sparse axpy */ 1498e745005bSHong Zhang nexta = 0; 1499d6ab1dc0SHong Zhang valtmp = pA[j]; 1500e745005bSHong Zhang for (k=0; nexta<anz; k++) { 1501e745005bSHong Zhang if (cj[k] == adj[nexta]) { 1502e745005bSHong Zhang ca[k] += valtmp*ada[nexta]; 1503e745005bSHong Zhang nexta++; 1504d6ab1dc0SHong Zhang } 1505e745005bSHong Zhang } 1506e745005bSHong Zhang ierr = PetscLogFlops(2.0*anz);CHKERRQ(ierr); 1507d6ab1dc0SHong Zhang } 1508d6ab1dc0SHong Zhang 1509d6ab1dc0SHong Zhang /* put the value into Cd (diagonal part) */ 1510d6ab1dc0SHong Zhang pnz = pd->i[i+1] - pd->i[i]; 1511d6ab1dc0SHong Zhang pdJ = pd->j + pd->i[i]; 1512d6ab1dc0SHong Zhang pA = pd->a + pd->i[i]; 1513d6ab1dc0SHong Zhang for (j=0; j<pnz; j++) { 1514d6ab1dc0SHong Zhang row = pdJ[j]; 1515d6ab1dc0SHong Zhang cj = bj + bi[row]; 1516d6ab1dc0SHong Zhang ca = ba + bi[row]; 1517e745005bSHong Zhang /* perform sparse axpy */ 1518e745005bSHong Zhang nexta = 0; 1519d6ab1dc0SHong Zhang valtmp = pA[j]; 1520e745005bSHong Zhang for (k=0; nexta<anz; k++) { 1521e745005bSHong Zhang if (cj[k] == adj[nexta]) { 1522e745005bSHong Zhang ca[k] += valtmp*ada[nexta]; 1523e745005bSHong Zhang nexta++; 1524d6ab1dc0SHong Zhang } 1525e745005bSHong Zhang } 1526e745005bSHong Zhang ierr = PetscLogFlops(2.0*anz);CHKERRQ(ierr); 1527d6ab1dc0SHong Zhang } 1528d6ab1dc0SHong Zhang } 1529d6ab1dc0SHong Zhang 1530d6ab1dc0SHong Zhang /* 3) send and recv matrix values coa */ 1531d6ab1dc0SHong Zhang /*------------------------------------*/ 1532d6ab1dc0SHong Zhang buf_ri = merge->buf_ri; 1533d6ab1dc0SHong Zhang buf_rj = merge->buf_rj; 1534d6ab1dc0SHong Zhang len_s = merge->len_s; 1535d6ab1dc0SHong Zhang ierr = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr); 1536d6ab1dc0SHong Zhang ierr = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr); 1537d6ab1dc0SHong Zhang 1538dcca6d9dSJed Brown ierr = PetscMalloc2(merge->nsend+1,&s_waits,size,&status);CHKERRQ(ierr); 1539d6ab1dc0SHong Zhang for (proc=0,k=0; proc<size; proc++) { 1540d6ab1dc0SHong Zhang if (!len_s[proc]) continue; 1541d6ab1dc0SHong Zhang i = merge->owners_co[proc]; 1542d6ab1dc0SHong Zhang ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr); 1543d6ab1dc0SHong Zhang k++; 1544d6ab1dc0SHong Zhang } 1545d6ab1dc0SHong Zhang if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);} 1546d6ab1dc0SHong Zhang if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);} 1547d6ab1dc0SHong Zhang 1548d6ab1dc0SHong Zhang ierr = PetscFree2(s_waits,status);CHKERRQ(ierr); 1549d6ab1dc0SHong Zhang ierr = PetscFree(r_waits);CHKERRQ(ierr); 1550d6ab1dc0SHong Zhang ierr = PetscFree(coa);CHKERRQ(ierr); 1551d6ab1dc0SHong Zhang 1552d6ab1dc0SHong Zhang /* 4) insert local Cseq and received values into Cmpi */ 1553d6ab1dc0SHong Zhang /*----------------------------------------------------*/ 1554dcca6d9dSJed Brown ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr); 1555d6ab1dc0SHong Zhang for (k=0; k<merge->nrecv; k++) { 1556e745005bSHong Zhang buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 1557d6ab1dc0SHong Zhang nrows = *(buf_ri_k[k]); 1558d6ab1dc0SHong Zhang nextrow[k] = buf_ri_k[k]+1; /* next row number of k-th recved i-structure */ 1559d6ab1dc0SHong Zhang nextci[k] = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure */ 1560d6ab1dc0SHong Zhang } 1561d6ab1dc0SHong Zhang 1562d6ab1dc0SHong Zhang for (i=0; i<cm; i++) { 1563d6ab1dc0SHong Zhang row = owners[rank] + i; /* global row index of C_seq */ 1564d6ab1dc0SHong Zhang bj_i = bj + bi[i]; /* col indices of the i-th row of C */ 1565d6ab1dc0SHong Zhang ba_i = ba + bi[i]; 1566d6ab1dc0SHong Zhang bnz = bi[i+1] - bi[i]; 1567d6ab1dc0SHong Zhang /* add received vals into ba */ 1568d6ab1dc0SHong Zhang for (k=0; k<merge->nrecv; k++) { /* k-th received message */ 1569d6ab1dc0SHong Zhang /* i-th row */ 1570d6ab1dc0SHong Zhang if (i == *nextrow[k]) { 1571d6ab1dc0SHong Zhang cnz = *(nextci[k]+1) - *nextci[k]; 1572d6ab1dc0SHong Zhang cj = buf_rj[k] + *(nextci[k]); 1573d6ab1dc0SHong Zhang ca = abuf_r[k] + *(nextci[k]); 1574d6ab1dc0SHong Zhang nextcj = 0; 1575d6ab1dc0SHong Zhang for (j=0; nextcj<cnz; j++) { 1576d6ab1dc0SHong Zhang if (bj_i[j] == cj[nextcj]) { /* bcol == ccol */ 1577d6ab1dc0SHong Zhang ba_i[j] += ca[nextcj++]; 1578d6ab1dc0SHong Zhang } 1579d6ab1dc0SHong Zhang } 1580d6ab1dc0SHong Zhang nextrow[k]++; nextci[k]++; 1581d6ab1dc0SHong Zhang ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr); 1582d6ab1dc0SHong Zhang } 1583d6ab1dc0SHong Zhang } 1584d6ab1dc0SHong Zhang ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr); 1585d6ab1dc0SHong Zhang } 1586d6ab1dc0SHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1587d6ab1dc0SHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1588d6ab1dc0SHong Zhang 1589d6ab1dc0SHong Zhang ierr = PetscFree(ba);CHKERRQ(ierr); 1590d6ab1dc0SHong Zhang ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr); 1591d6ab1dc0SHong Zhang ierr = PetscFree(abuf_r);CHKERRQ(ierr); 1592d6ab1dc0SHong Zhang ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr); 1593d6ab1dc0SHong Zhang PetscFunctionReturn(0); 1594d6ab1dc0SHong Zhang } 1595d6ab1dc0SHong Zhang 1596c9ba551fSBarry Smith PetscErrorCode MatDuplicate_MPIAIJ_MatPtAP(Mat, MatDuplicateOption,Mat*); 15972bbb1c24SHong Zhang /* This routine is modified from MatPtAPSymbolic_MPIAIJ_MPIAIJ(); 15982bbb1c24SHong Zhang differ from MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable in using LLCondensedCreate_Scalable() */ 1599d6ab1dc0SHong Zhang #undef __FUNCT__ 16006da69ca6SHong Zhang #define __FUNCT__ "MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ" 16016da69ca6SHong Zhang PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(Mat P,Mat A,PetscReal fill,Mat *C) 1602d6ab1dc0SHong Zhang { 1603d6ab1dc0SHong Zhang PetscErrorCode ierr; 1604d6ab1dc0SHong Zhang Mat Cmpi,A_loc,POt,PDt; 1605d6ab1dc0SHong Zhang Mat_PtAPMPI *ptap; 16060298fd71SBarry Smith PetscFreeSpaceList free_space=NULL,current_space=NULL; 1607d6ab1dc0SHong Zhang Mat_MPIAIJ *p =(Mat_MPIAIJ*)P->data,*c; 1608d6ab1dc0SHong Zhang PetscInt *pdti,*pdtj,*poti,*potj,*ptJ; 1609d6ab1dc0SHong Zhang PetscInt nnz; 1610d6ab1dc0SHong Zhang PetscInt *lnk,*owners_co,*coi,*coj,i,k,pnz,row; 1611d6ab1dc0SHong Zhang PetscInt am =A->rmap->n,pn=P->cmap->n; 1612ce94432eSBarry Smith MPI_Comm comm; 1613d6ab1dc0SHong Zhang PetscMPIInt size,rank,tagi,tagj,*len_si,*len_s,*len_ri; 1614d6ab1dc0SHong Zhang PetscInt **buf_rj,**buf_ri,**buf_ri_k; 1615d6ab1dc0SHong Zhang PetscInt len,proc,*dnz,*onz,*owners; 1616d6ab1dc0SHong Zhang PetscInt nzi,*bi,*bj; 1617d6ab1dc0SHong Zhang PetscInt nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci; 1618d6ab1dc0SHong Zhang MPI_Request *swaits,*rwaits; 1619d6ab1dc0SHong Zhang MPI_Status *sstatus,rstatus; 1620d6ab1dc0SHong Zhang Mat_Merge_SeqsToMPI *merge; 1621d6ab1dc0SHong Zhang PetscInt *ai,*aj,*Jptr,anz,*prmap=p->garray,pon,nspacedouble=0,j; 1622d6ab1dc0SHong Zhang PetscReal afill =1.0,afill_tmp; 1623c36aecf5SHong Zhang PetscInt rstart = P->cmap->rstart,rmax,aN=A->cmap->N,Crmax; 1624d6ab1dc0SHong Zhang PetscScalar *vals; 1625d6ab1dc0SHong Zhang Mat_SeqAIJ *a_loc, *pdt,*pot; 1626d6ab1dc0SHong Zhang 1627d6ab1dc0SHong Zhang PetscFunctionBegin; 1628ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 1629d6ab1dc0SHong Zhang /* check if matrix local sizes are compatible */ 1630d6ab1dc0SHong Zhang if (A->rmap->rstart != P->rmap->rstart || A->rmap->rend != P->rmap->rend) { 1631d6ab1dc0SHong 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); 1632d6ab1dc0SHong Zhang } 1633d6ab1dc0SHong Zhang 1634d6ab1dc0SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1635d6ab1dc0SHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 1636d6ab1dc0SHong Zhang 1637d6ab1dc0SHong Zhang /* create struct Mat_PtAPMPI and attached it to C later */ 1638b00a9115SJed Brown ierr = PetscNew(&ptap);CHKERRQ(ierr); 1639d6ab1dc0SHong Zhang 1640d6ab1dc0SHong Zhang /* get A_loc by taking all local rows of A */ 1641d6ab1dc0SHong Zhang ierr = MatMPIAIJGetLocalMat(A,MAT_INITIAL_MATRIX,&A_loc);CHKERRQ(ierr); 16422205254eSKarl Rupp 1643d6ab1dc0SHong Zhang ptap->A_loc = A_loc; 1644d6ab1dc0SHong Zhang a_loc = (Mat_SeqAIJ*)(A_loc)->data; 1645d6ab1dc0SHong Zhang ai = a_loc->i; 1646d6ab1dc0SHong Zhang aj = a_loc->j; 1647d6ab1dc0SHong Zhang 1648d6ab1dc0SHong Zhang /* determine symbolic Co=(p->B)^T*A - send to others */ 1649d6ab1dc0SHong Zhang /*----------------------------------------------------*/ 1650d6ab1dc0SHong Zhang ierr = MatTransposeSymbolic_SeqAIJ(p->A,&PDt);CHKERRQ(ierr); 1651d6ab1dc0SHong Zhang pdt = (Mat_SeqAIJ*)PDt->data; 1652d6ab1dc0SHong Zhang pdti = pdt->i; pdtj = pdt->j; 1653d6ab1dc0SHong Zhang 1654d6ab1dc0SHong Zhang ierr = MatTransposeSymbolic_SeqAIJ(p->B,&POt);CHKERRQ(ierr); 1655d6ab1dc0SHong Zhang pot = (Mat_SeqAIJ*)POt->data; 1656d6ab1dc0SHong Zhang poti = pot->i; potj = pot->j; 1657d6ab1dc0SHong Zhang 1658d6ab1dc0SHong Zhang /* then, compute symbolic Co = (p->B)^T*A */ 1659d6ab1dc0SHong Zhang pon = (p->B)->cmap->n; /* total num of rows to be sent to other processors 1660d6ab1dc0SHong Zhang >= (num of nonzero rows of C_seq) - pn */ 1661854ce69bSBarry Smith ierr = PetscMalloc1(pon+1,&coi);CHKERRQ(ierr); 1662d6ab1dc0SHong Zhang coi[0] = 0; 1663d6ab1dc0SHong Zhang 1664d6ab1dc0SHong Zhang /* set initial free space to be fill*(nnz(p->B) + nnz(A)) */ 1665d6ab1dc0SHong Zhang nnz = fill*(poti[pon] + ai[am]); 1666a2ea699eSBarry Smith ierr = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr); 1667d6ab1dc0SHong Zhang current_space = free_space; 166819f0e57aSHong Zhang 1669d6ab1dc0SHong Zhang /* create and initialize a linked list */ 1670d6ab1dc0SHong Zhang i = PetscMax(pdt->rmax,pot->rmax); 1671c36aecf5SHong Zhang Crmax = i*a_loc->rmax*size; /* non-scalable! */ 1672d6ab1dc0SHong Zhang if (!Crmax || Crmax > aN) Crmax = aN; 1673d6ab1dc0SHong Zhang ierr = PetscLLCondensedCreate_Scalable(Crmax,&lnk);CHKERRQ(ierr); 1674d6ab1dc0SHong Zhang 1675d6ab1dc0SHong Zhang for (i=0; i<pon; i++) { 1676d6ab1dc0SHong Zhang pnz = poti[i+1] - poti[i]; 1677d6ab1dc0SHong Zhang ptJ = potj + poti[i]; 1678d6ab1dc0SHong Zhang for (j=0; j<pnz; j++) { 1679d6ab1dc0SHong Zhang row = ptJ[j]; /* row of A_loc == col of Pot */ 1680d6ab1dc0SHong Zhang anz = ai[row+1] - ai[row]; 1681d6ab1dc0SHong Zhang Jptr = aj + ai[row]; 1682d6ab1dc0SHong Zhang /* add non-zero cols of AP into the sorted linked list lnk */ 1683d6ab1dc0SHong Zhang ierr = PetscLLCondensedAddSorted_Scalable(anz,Jptr,lnk);CHKERRQ(ierr); 1684d6ab1dc0SHong Zhang } 1685d6ab1dc0SHong Zhang nnz = lnk[0]; 1686d6ab1dc0SHong Zhang 1687d6ab1dc0SHong Zhang /* If free space is not available, double the total space in the list */ 1688d6ab1dc0SHong Zhang if (current_space->local_remaining<nnz) { 1689d6ab1dc0SHong Zhang ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 1690d6ab1dc0SHong Zhang nspacedouble++; 1691d6ab1dc0SHong Zhang } 1692d6ab1dc0SHong Zhang 1693d6ab1dc0SHong Zhang /* Copy data into free space, and zero out denserows */ 1694d6ab1dc0SHong Zhang ierr = PetscLLCondensedClean_Scalable(nnz,current_space->array,lnk);CHKERRQ(ierr); 16952205254eSKarl Rupp 1696d6ab1dc0SHong Zhang current_space->array += nnz; 1697d6ab1dc0SHong Zhang current_space->local_used += nnz; 1698d6ab1dc0SHong Zhang current_space->local_remaining -= nnz; 16992205254eSKarl Rupp 1700d6ab1dc0SHong Zhang coi[i+1] = coi[i] + nnz; 1701d6ab1dc0SHong Zhang } 1702d6ab1dc0SHong Zhang 1703854ce69bSBarry Smith ierr = PetscMalloc1(coi[pon]+1,&coj);CHKERRQ(ierr); 1704d6ab1dc0SHong Zhang ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr); 17052205254eSKarl Rupp 1706118727c9SMark F. Adams afill_tmp = (PetscReal)coi[pon]/(poti[pon] + ai[am]+1); 1707d6ab1dc0SHong Zhang if (afill_tmp > afill) afill = afill_tmp; 1708d6ab1dc0SHong Zhang 1709d6ab1dc0SHong Zhang /* send j-array (coj) of Co to other processors */ 1710d6ab1dc0SHong Zhang /*----------------------------------------------*/ 1711d6ab1dc0SHong Zhang /* determine row ownership */ 1712b00a9115SJed Brown ierr = PetscNew(&merge);CHKERRQ(ierr); 1713d6ab1dc0SHong Zhang ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr); 17142205254eSKarl Rupp 1715d6ab1dc0SHong Zhang merge->rowmap->n = pn; 1716d6ab1dc0SHong Zhang merge->rowmap->bs = 1; 17172205254eSKarl Rupp 1718d6ab1dc0SHong Zhang ierr = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr); 1719d6ab1dc0SHong Zhang owners = merge->rowmap->range; 1720d6ab1dc0SHong Zhang 1721d6ab1dc0SHong Zhang /* determine the number of messages to send, their lengths */ 17221795a4d1SJed Brown ierr = PetscCalloc1(size,&len_si);CHKERRQ(ierr); 1723785e854fSJed Brown ierr = PetscMalloc1(size,&merge->len_s);CHKERRQ(ierr); 17242205254eSKarl Rupp 1725d6ab1dc0SHong Zhang len_s = merge->len_s; 1726d6ab1dc0SHong Zhang merge->nsend = 0; 1727d6ab1dc0SHong Zhang 1728854ce69bSBarry Smith ierr = PetscMalloc1(size+2,&owners_co);CHKERRQ(ierr); 1729d6ab1dc0SHong Zhang ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 1730d6ab1dc0SHong Zhang 1731d6ab1dc0SHong Zhang proc = 0; 1732d6ab1dc0SHong Zhang for (i=0; i<pon; i++) { 1733d6ab1dc0SHong Zhang while (prmap[i] >= owners[proc+1]) proc++; 1734d6ab1dc0SHong Zhang len_si[proc]++; /* num of rows in Co to be sent to [proc] */ 1735d6ab1dc0SHong Zhang len_s[proc] += coi[i+1] - coi[i]; 1736d6ab1dc0SHong Zhang } 1737d6ab1dc0SHong Zhang 1738d6ab1dc0SHong Zhang len = 0; /* max length of buf_si[] */ 1739d6ab1dc0SHong Zhang owners_co[0] = 0; 1740d6ab1dc0SHong Zhang for (proc=0; proc<size; proc++) { 1741d6ab1dc0SHong Zhang owners_co[proc+1] = owners_co[proc] + len_si[proc]; 1742d6ab1dc0SHong Zhang if (len_si[proc]) { 1743d6ab1dc0SHong Zhang merge->nsend++; 1744d6ab1dc0SHong Zhang len_si[proc] = 2*(len_si[proc] + 1); 1745d6ab1dc0SHong Zhang len += len_si[proc]; 1746d6ab1dc0SHong Zhang } 1747d6ab1dc0SHong Zhang } 1748d6ab1dc0SHong Zhang 1749d6ab1dc0SHong Zhang /* determine the number and length of messages to receive for coi and coj */ 17500298fd71SBarry Smith ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&merge->nrecv);CHKERRQ(ierr); 1751d6ab1dc0SHong Zhang ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr); 1752d6ab1dc0SHong Zhang 1753d6ab1dc0SHong Zhang /* post the Irecv and Isend of coj */ 1754d6ab1dc0SHong Zhang ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr); 1755d6ab1dc0SHong Zhang ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr); 1756854ce69bSBarry Smith ierr = PetscMalloc1(merge->nsend+1,&swaits);CHKERRQ(ierr); 1757d6ab1dc0SHong Zhang for (proc=0, k=0; proc<size; proc++) { 1758d6ab1dc0SHong Zhang if (!len_s[proc]) continue; 1759d6ab1dc0SHong Zhang i = owners_co[proc]; 1760d6ab1dc0SHong Zhang ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr); 1761d6ab1dc0SHong Zhang k++; 1762d6ab1dc0SHong Zhang } 1763d6ab1dc0SHong Zhang 1764d6ab1dc0SHong Zhang /* receives and sends of coj are complete */ 1765785e854fSJed Brown ierr = PetscMalloc1(size,&sstatus);CHKERRQ(ierr); 1766d6ab1dc0SHong Zhang for (i=0; i<merge->nrecv; i++) { 1767c280ed6aSJed Brown PetscMPIInt icompleted; 1768c280ed6aSJed Brown ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 1769d6ab1dc0SHong Zhang } 1770d6ab1dc0SHong Zhang ierr = PetscFree(rwaits);CHKERRQ(ierr); 1771d6ab1dc0SHong Zhang if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);} 1772d6ab1dc0SHong Zhang 1773d6ab1dc0SHong Zhang /* send and recv coi */ 1774d6ab1dc0SHong Zhang /*-------------------*/ 1775d6ab1dc0SHong Zhang ierr = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr); 1776d6ab1dc0SHong Zhang ierr = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr); 1777854ce69bSBarry Smith ierr = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr); 1778d6ab1dc0SHong Zhang buf_si = buf_s; /* points to the beginning of k-th msg to be sent */ 1779d6ab1dc0SHong Zhang for (proc=0,k=0; proc<size; proc++) { 1780d6ab1dc0SHong Zhang if (!len_s[proc]) continue; 1781d6ab1dc0SHong Zhang /* form outgoing message for i-structure: 1782d6ab1dc0SHong Zhang buf_si[0]: nrows to be sent 1783d6ab1dc0SHong Zhang [1:nrows]: row index (global) 1784d6ab1dc0SHong Zhang [nrows+1:2*nrows+1]: i-structure index 1785d6ab1dc0SHong Zhang */ 1786d6ab1dc0SHong Zhang /*-------------------------------------------*/ 1787d6ab1dc0SHong Zhang nrows = len_si[proc]/2 - 1; 1788d6ab1dc0SHong Zhang buf_si_i = buf_si + nrows+1; 1789d6ab1dc0SHong Zhang buf_si[0] = nrows; 1790d6ab1dc0SHong Zhang buf_si_i[0] = 0; 1791d6ab1dc0SHong Zhang nrows = 0; 1792d6ab1dc0SHong Zhang for (i=owners_co[proc]; i<owners_co[proc+1]; i++) { 1793d6ab1dc0SHong Zhang nzi = coi[i+1] - coi[i]; 1794d6ab1dc0SHong Zhang buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */ 1795d6ab1dc0SHong Zhang buf_si[nrows+1] = prmap[i] -owners[proc]; /* local row index */ 1796d6ab1dc0SHong Zhang nrows++; 1797d6ab1dc0SHong Zhang } 1798d6ab1dc0SHong Zhang ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr); 1799d6ab1dc0SHong Zhang k++; 1800d6ab1dc0SHong Zhang buf_si += len_si[proc]; 1801d6ab1dc0SHong Zhang } 1802d6ab1dc0SHong Zhang i = merge->nrecv; 1803d6ab1dc0SHong Zhang while (i--) { 1804c280ed6aSJed Brown PetscMPIInt icompleted; 1805c280ed6aSJed Brown ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 1806d6ab1dc0SHong Zhang } 1807d6ab1dc0SHong Zhang ierr = PetscFree(rwaits);CHKERRQ(ierr); 1808d6ab1dc0SHong Zhang if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);} 1809d6ab1dc0SHong Zhang ierr = PetscFree(len_si);CHKERRQ(ierr); 1810d6ab1dc0SHong Zhang ierr = PetscFree(len_ri);CHKERRQ(ierr); 1811d6ab1dc0SHong Zhang ierr = PetscFree(swaits);CHKERRQ(ierr); 1812d6ab1dc0SHong Zhang ierr = PetscFree(sstatus);CHKERRQ(ierr); 1813d6ab1dc0SHong Zhang ierr = PetscFree(buf_s);CHKERRQ(ierr); 1814d6ab1dc0SHong Zhang 1815d6ab1dc0SHong Zhang /* compute the local portion of C (mpi mat) */ 1816d6ab1dc0SHong Zhang /*------------------------------------------*/ 1817d6ab1dc0SHong Zhang /* allocate bi array and free space for accumulating nonzero column info */ 1818854ce69bSBarry Smith ierr = PetscMalloc1(pn+1,&bi);CHKERRQ(ierr); 1819d6ab1dc0SHong Zhang bi[0] = 0; 1820d6ab1dc0SHong Zhang 1821d6ab1dc0SHong Zhang /* set initial free space to be fill*(nnz(P) + nnz(AP)) */ 1822d6ab1dc0SHong Zhang nnz = fill*(pdti[pn] + poti[pon] + ai[am]); 1823a2ea699eSBarry Smith ierr = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr); 1824d6ab1dc0SHong Zhang current_space = free_space; 1825d6ab1dc0SHong Zhang 1826dcca6d9dSJed Brown ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr); 1827d6ab1dc0SHong Zhang for (k=0; k<merge->nrecv; k++) { 1828d6ab1dc0SHong Zhang buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 1829d6ab1dc0SHong Zhang nrows = *buf_ri_k[k]; 1830d6ab1dc0SHong Zhang nextrow[k] = buf_ri_k[k] + 1; /* next row number of k-th recved i-structure */ 18312205254eSKarl Rupp nextci[k] = buf_ri_k[k] + (nrows + 1); /* points to the next i-structure of k-th recieved i-structure */ 1832d6ab1dc0SHong Zhang } 1833d6ab1dc0SHong Zhang 1834d6ab1dc0SHong Zhang ierr = MatPreallocateInitialize(comm,pn,A->cmap->n,dnz,onz);CHKERRQ(ierr); 1835d6ab1dc0SHong Zhang rmax = 0; 1836d6ab1dc0SHong Zhang for (i=0; i<pn; i++) { 1837d6ab1dc0SHong Zhang /* add pdt[i,:]*AP into lnk */ 1838d6ab1dc0SHong Zhang pnz = pdti[i+1] - pdti[i]; 1839d6ab1dc0SHong Zhang ptJ = pdtj + pdti[i]; 1840d6ab1dc0SHong Zhang for (j=0; j<pnz; j++) { 1841d6ab1dc0SHong Zhang row = ptJ[j]; /* row of AP == col of Pt */ 1842d6ab1dc0SHong Zhang anz = ai[row+1] - ai[row]; 1843d6ab1dc0SHong Zhang Jptr = aj + ai[row]; 1844d6ab1dc0SHong Zhang /* add non-zero cols of AP into the sorted linked list lnk */ 1845d6ab1dc0SHong Zhang ierr = PetscLLCondensedAddSorted_Scalable(anz,Jptr,lnk);CHKERRQ(ierr); 1846d6ab1dc0SHong Zhang } 1847d6ab1dc0SHong Zhang 1848d6ab1dc0SHong Zhang /* add received col data into lnk */ 1849d6ab1dc0SHong Zhang for (k=0; k<merge->nrecv; k++) { /* k-th received message */ 1850d6ab1dc0SHong Zhang if (i == *nextrow[k]) { /* i-th row */ 1851d6ab1dc0SHong Zhang nzi = *(nextci[k]+1) - *nextci[k]; 1852d6ab1dc0SHong Zhang Jptr = buf_rj[k] + *nextci[k]; 1853d6ab1dc0SHong Zhang ierr = PetscLLCondensedAddSorted_Scalable(nzi,Jptr,lnk);CHKERRQ(ierr); 1854d6ab1dc0SHong Zhang nextrow[k]++; nextci[k]++; 1855d6ab1dc0SHong Zhang } 1856d6ab1dc0SHong Zhang } 1857d6ab1dc0SHong Zhang nnz = lnk[0]; 1858d6ab1dc0SHong Zhang 1859d6ab1dc0SHong Zhang /* if free space is not available, make more free space */ 1860d6ab1dc0SHong Zhang if (current_space->local_remaining<nnz) { 1861d6ab1dc0SHong Zhang ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 1862d6ab1dc0SHong Zhang nspacedouble++; 1863d6ab1dc0SHong Zhang } 1864d6ab1dc0SHong Zhang /* copy data into free space, then initialize lnk */ 1865d6ab1dc0SHong Zhang ierr = PetscLLCondensedClean_Scalable(nnz,current_space->array,lnk);CHKERRQ(ierr); 1866d6ab1dc0SHong Zhang ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr); 18672205254eSKarl Rupp 1868d6ab1dc0SHong Zhang current_space->array += nnz; 1869d6ab1dc0SHong Zhang current_space->local_used += nnz; 1870d6ab1dc0SHong Zhang current_space->local_remaining -= nnz; 18712205254eSKarl Rupp 1872d6ab1dc0SHong Zhang bi[i+1] = bi[i] + nnz; 1873d6ab1dc0SHong Zhang if (nnz > rmax) rmax = nnz; 1874d6ab1dc0SHong Zhang } 1875f671be37SHong Zhang ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr); 1876d6ab1dc0SHong Zhang 1877854ce69bSBarry Smith ierr = PetscMalloc1(bi[pn]+1,&bj);CHKERRQ(ierr); 1878d6ab1dc0SHong Zhang ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr); 1879118727c9SMark F. Adams afill_tmp = (PetscReal)bi[pn]/(pdti[pn] + poti[pon] + ai[am]+1); 1880d6ab1dc0SHong Zhang if (afill_tmp > afill) afill = afill_tmp; 1881d6ab1dc0SHong Zhang ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr); 1882d6ab1dc0SHong Zhang ierr = MatDestroy(&POt);CHKERRQ(ierr); 1883d6ab1dc0SHong Zhang ierr = MatDestroy(&PDt);CHKERRQ(ierr); 1884d6ab1dc0SHong Zhang 1885d6ab1dc0SHong Zhang /* create symbolic parallel matrix Cmpi - why cannot be assembled in Numeric part */ 1886d6ab1dc0SHong Zhang /*----------------------------------------------------------------------------------*/ 18871795a4d1SJed Brown ierr = PetscCalloc1(rmax+1,&vals);CHKERRQ(ierr); 1888d6ab1dc0SHong Zhang 1889d6ab1dc0SHong Zhang ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr); 1890d6ab1dc0SHong Zhang ierr = MatSetSizes(Cmpi,pn,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 189133d57670SJed Brown ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(A->cmap->bs));CHKERRQ(ierr); 1892d6ab1dc0SHong Zhang ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr); 1893d6ab1dc0SHong Zhang ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr); 1894d6ab1dc0SHong Zhang ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 1895d6ab1dc0SHong Zhang ierr = MatSetBlockSize(Cmpi,1);CHKERRQ(ierr); 1896d6ab1dc0SHong Zhang for (i=0; i<pn; i++) { 1897d6ab1dc0SHong Zhang row = i + rstart; 1898d6ab1dc0SHong Zhang nnz = bi[i+1] - bi[i]; 1899d6ab1dc0SHong Zhang Jptr = bj + bi[i]; 1900d6ab1dc0SHong Zhang ierr = MatSetValues(Cmpi,1,&row,nnz,Jptr,vals,INSERT_VALUES);CHKERRQ(ierr); 1901d6ab1dc0SHong Zhang } 1902d6ab1dc0SHong Zhang ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1903d6ab1dc0SHong Zhang ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1904d6ab1dc0SHong Zhang ierr = PetscFree(vals);CHKERRQ(ierr); 1905d6ab1dc0SHong Zhang 1906d6ab1dc0SHong Zhang merge->bi = bi; 1907d6ab1dc0SHong Zhang merge->bj = bj; 1908d6ab1dc0SHong Zhang merge->coi = coi; 1909d6ab1dc0SHong Zhang merge->coj = coj; 1910d6ab1dc0SHong Zhang merge->buf_ri = buf_ri; 1911d6ab1dc0SHong Zhang merge->buf_rj = buf_rj; 1912d6ab1dc0SHong Zhang merge->owners_co = owners_co; 1913d6ab1dc0SHong Zhang merge->destroy = Cmpi->ops->destroy; 1914d6ab1dc0SHong Zhang merge->duplicate = Cmpi->ops->duplicate; 1915d6ab1dc0SHong Zhang 19166da69ca6SHong Zhang Cmpi->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ; 1917d6ab1dc0SHong Zhang Cmpi->ops->destroy = MatDestroy_MPIAIJ_PtAP; 1918c9ba551fSBarry Smith Cmpi->ops->duplicate = MatDuplicate_MPIAIJ_MatPtAP; 1919d6ab1dc0SHong Zhang 1920d6ab1dc0SHong Zhang /* attach the supporting struct to Cmpi for reuse */ 1921d6ab1dc0SHong Zhang c = (Mat_MPIAIJ*)Cmpi->data; 19222205254eSKarl Rupp 1923d6ab1dc0SHong Zhang c->ptap = ptap; 19240298fd71SBarry Smith ptap->api = NULL; 19250298fd71SBarry Smith ptap->apj = NULL; 1926d6ab1dc0SHong Zhang ptap->merge = merge; 1927d6ab1dc0SHong Zhang ptap->rmax = rmax; 19280298fd71SBarry Smith ptap->apa = NULL; 1929f96d37fcSHong Zhang 1930f96d37fcSHong Zhang *C = Cmpi; 1931f96d37fcSHong Zhang #if defined(PETSC_USE_INFO) 1932f96d37fcSHong Zhang if (bi[pn] != 0) { 193357622a8eSBarry Smith ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr); 193457622a8eSBarry Smith ierr = PetscInfo1(Cmpi,"Use MatTransposeMatMult(A,B,MatReuse,%g,&C) for best performance.\n",(double)afill);CHKERRQ(ierr); 1935f96d37fcSHong Zhang } else { 1936f96d37fcSHong Zhang ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr); 1937f96d37fcSHong Zhang } 1938f96d37fcSHong Zhang #endif 1939f96d37fcSHong Zhang PetscFunctionReturn(0); 1940f96d37fcSHong Zhang } 1941