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"}; 190fc8cf34SHong Zhang PetscInt alg=0; /* 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; 90598bc09dSHong Zhang PetscInt *adi=ad->i,*adj,*aoi=ao->i,*aoj; 91598bc09dSHong Zhang PetscScalar *ada,*aoa,*cda=cd->a,*coa=co->a; 92598bc09dSHong Zhang Mat_SeqAIJ *p_loc,*p_oth; 93598bc09dSHong Zhang PetscInt *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pj; 94598bc09dSHong Zhang PetscScalar *pa_loc,*pa_oth,*pa,*apa,valtmp,*ca; 95598bc09dSHong Zhang PetscInt cm =C->rmap->n,anz,pnz; 96598bc09dSHong Zhang Mat_PtAPMPI *ptap=c->ptap; 9765a57a32SSean Farley PetscInt *api,*apj,*apJ,i,j,k,row; 9829825010SHong Zhang PetscInt cstart=C->cmap->rstart; 99598bc09dSHong Zhang PetscInt cdnz,conz,k0,k1; 1009467f45fSHong Zhang MPI_Comm comm; 1019467f45fSHong Zhang PetscMPIInt size; 102598bc09dSHong Zhang 103598bc09dSHong Zhang PetscFunctionBegin; 1049467f45fSHong Zhang ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 1059467f45fSHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1069467f45fSHong Zhang 107a1a4e84aSHong Zhang /* 1) get P_oth = ptap->P_oth and P_loc = ptap->P_loc */ 108598bc09dSHong Zhang /*-----------------------------------------------------*/ 109a1a4e84aSHong Zhang /* update numerical values of P_oth and P_loc */ 110b7f45c76SHong Zhang ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 111a1a4e84aSHong Zhang ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 112598bc09dSHong Zhang 113598bc09dSHong Zhang /* 2) compute numeric C_loc = A_loc*P = Ad*P_loc + Ao*P_oth */ 114598bc09dSHong Zhang /*----------------------------------------------------------*/ 115598bc09dSHong Zhang /* get data from symbolic products */ 116a1a4e84aSHong Zhang p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 117598bc09dSHong Zhang pi_loc=p_loc->i; pj_loc=p_loc->j; pa_loc=p_loc->a; 1189467f45fSHong Zhang if (size >1) { 1199467f45fSHong Zhang p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 120598bc09dSHong Zhang pi_oth=p_oth->i; pj_oth=p_oth->j; pa_oth=p_oth->a; 121968382fdSHong Zhang } else { 122968382fdSHong Zhang pi_oth=NULL; pj_oth=NULL; pa_oth=NULL; 1239467f45fSHong Zhang } 124598bc09dSHong Zhang 125598bc09dSHong Zhang /* get apa for storing dense row A[i,:]*P */ 126598bc09dSHong Zhang apa = ptap->apa; 127598bc09dSHong Zhang 12829825010SHong Zhang api = ptap->api; 12929825010SHong Zhang apj = ptap->apj; 130598bc09dSHong Zhang for (i=0; i<cm; i++) { 131598bc09dSHong Zhang /* diagonal portion of A */ 132598bc09dSHong Zhang anz = adi[i+1] - adi[i]; 133598bc09dSHong Zhang adj = ad->j + adi[i]; 134598bc09dSHong Zhang ada = ad->a + adi[i]; 135598bc09dSHong Zhang for (j=0; j<anz; j++) { 136598bc09dSHong Zhang row = adj[j]; 137598bc09dSHong Zhang pnz = pi_loc[row+1] - pi_loc[row]; 138598bc09dSHong Zhang pj = pj_loc + pi_loc[row]; 139598bc09dSHong Zhang pa = pa_loc + pi_loc[row]; 140598bc09dSHong Zhang 141598bc09dSHong Zhang /* perform dense axpy */ 142598bc09dSHong Zhang valtmp = ada[j]; 143598bc09dSHong Zhang for (k=0; k<pnz; k++) { 144598bc09dSHong Zhang apa[pj[k]] += valtmp*pa[k]; 145598bc09dSHong Zhang } 146598bc09dSHong Zhang ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr); 147598bc09dSHong Zhang } 148598bc09dSHong Zhang 149598bc09dSHong Zhang /* off-diagonal portion of A */ 150598bc09dSHong Zhang anz = aoi[i+1] - aoi[i]; 151598bc09dSHong Zhang aoj = ao->j + aoi[i]; 152598bc09dSHong Zhang aoa = ao->a + aoi[i]; 153598bc09dSHong Zhang for (j=0; j<anz; j++) { 154598bc09dSHong Zhang row = aoj[j]; 155598bc09dSHong Zhang pnz = pi_oth[row+1] - pi_oth[row]; 156598bc09dSHong Zhang pj = pj_oth + pi_oth[row]; 157598bc09dSHong Zhang pa = pa_oth + pi_oth[row]; 158598bc09dSHong Zhang 159598bc09dSHong Zhang /* perform dense axpy */ 160598bc09dSHong Zhang valtmp = aoa[j]; 161598bc09dSHong Zhang for (k=0; k<pnz; k++) { 162598bc09dSHong Zhang apa[pj[k]] += valtmp*pa[k]; 163598bc09dSHong Zhang } 164598bc09dSHong Zhang ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr); 165598bc09dSHong Zhang } 166598bc09dSHong Zhang 167598bc09dSHong Zhang /* set values in C */ 168598bc09dSHong Zhang apJ = apj + api[i]; 169598bc09dSHong Zhang cdnz = cd->i[i+1] - cd->i[i]; 170598bc09dSHong Zhang conz = co->i[i+1] - co->i[i]; 171598bc09dSHong Zhang 172598bc09dSHong Zhang /* 1st off-diagoanl part of C */ 173598bc09dSHong Zhang ca = coa + co->i[i]; 174598bc09dSHong Zhang k = 0; 175598bc09dSHong Zhang for (k0=0; k0<conz; k0++) { 176598bc09dSHong Zhang if (apJ[k] >= cstart) break; 177598bc09dSHong Zhang ca[k0] = apa[apJ[k]]; 178598bc09dSHong Zhang apa[apJ[k]] = 0.0; 179598bc09dSHong Zhang k++; 180598bc09dSHong Zhang } 181598bc09dSHong Zhang 182598bc09dSHong Zhang /* diagonal part of C */ 183598bc09dSHong Zhang ca = cda + cd->i[i]; 184598bc09dSHong Zhang for (k1=0; k1<cdnz; k1++) { 185598bc09dSHong Zhang ca[k1] = apa[apJ[k]]; 186598bc09dSHong Zhang apa[apJ[k]] = 0.0; 187598bc09dSHong Zhang k++; 188598bc09dSHong Zhang } 189598bc09dSHong Zhang 190598bc09dSHong Zhang /* 2nd off-diagoanl part of C */ 191598bc09dSHong Zhang ca = coa + co->i[i]; 192598bc09dSHong Zhang for (; k0<conz; k0++) { 193598bc09dSHong Zhang ca[k0] = apa[apJ[k]]; 194598bc09dSHong Zhang apa[apJ[k]] = 0.0; 195598bc09dSHong Zhang k++; 196598bc09dSHong Zhang } 197598bc09dSHong Zhang } 198598bc09dSHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 199598bc09dSHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 200598bc09dSHong Zhang PetscFunctionReturn(0); 201598bc09dSHong Zhang } 202598bc09dSHong Zhang 203598bc09dSHong Zhang #undef __FUNCT__ 2040fc8cf34SHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable" 2050fc8cf34SHong Zhang PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat A,Mat P,PetscReal fill,Mat *C) 206598bc09dSHong Zhang { 207598bc09dSHong Zhang PetscErrorCode ierr; 208ce94432eSBarry Smith MPI_Comm comm; 2099467f45fSHong Zhang PetscMPIInt size; 210bfcd1a73SHong Zhang Mat Cmpi; 211598bc09dSHong Zhang Mat_PtAPMPI *ptap; 2120298fd71SBarry Smith PetscFreeSpaceList free_space=NULL,current_space=NULL; 2134ae0847bSHong Zhang Mat_MPIAIJ *a =(Mat_MPIAIJ*)A->data,*c; 214bfcd1a73SHong Zhang Mat_SeqAIJ *ad =(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_loc,*p_oth; 2154ae0847bSHong Zhang PetscInt *pi_loc,*pj_loc,*pi_oth,*pj_oth,*dnz,*onz; 2164ae0847bSHong Zhang PetscInt *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,rstart=A->rmap->rstart; 217d14fa243SHong Zhang PetscInt *lnk,i,pnz,row,*api,*apj,*Jptr,apnz,nspacedouble=0,j,nzi; 218bfcd1a73SHong Zhang PetscInt am=A->rmap->n,pN=P->cmap->N,pn=P->cmap->n,pm=P->rmap->n; 219598bc09dSHong Zhang PetscBT lnkbt; 220598bc09dSHong Zhang PetscScalar *apa; 221bfcd1a73SHong Zhang PetscReal afill; 222f84c536bSHong Zhang PetscInt nlnk_max,armax,prmax; 223598bc09dSHong Zhang 224598bc09dSHong Zhang PetscFunctionBegin; 225ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 2269467f45fSHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 2279467f45fSHong Zhang 228a1a4e84aSHong Zhang if (A->cmap->rstart != P->rmap->rstart || A->cmap->rend != P->rmap->rend) { 229cf1a0672SHong 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); 230a1a4e84aSHong Zhang } 231152983bfSHong Zhang 232598bc09dSHong Zhang /* create struct Mat_PtAPMPI and attached it to C later */ 233b00a9115SJed Brown ierr = PetscNew(&ptap);CHKERRQ(ierr); 234598bc09dSHong Zhang 235598bc09dSHong Zhang /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */ 236b7f45c76SHong Zhang ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 23719f0e57aSHong Zhang 238598bc09dSHong Zhang /* get P_loc by taking all local rows of P */ 239a1a4e84aSHong Zhang ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 240598bc09dSHong Zhang 241a1a4e84aSHong Zhang p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 242598bc09dSHong Zhang pi_loc = p_loc->i; pj_loc = p_loc->j; 2439467f45fSHong Zhang if (size > 1) { 2449467f45fSHong Zhang p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 245598bc09dSHong Zhang pi_oth = p_oth->i; pj_oth = p_oth->j; 24620e3dc75SHong Zhang } else { 24720e3dc75SHong Zhang p_oth = NULL; 24820e3dc75SHong Zhang pi_oth = NULL; pj_oth = NULL; 2499467f45fSHong Zhang } 250598bc09dSHong Zhang 251598bc09dSHong Zhang /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */ 252598bc09dSHong Zhang /*-------------------------------------------------------------------*/ 253854ce69bSBarry Smith ierr = PetscMalloc1(am+2,&api);CHKERRQ(ierr); 254a1a4e84aSHong Zhang ptap->api = api; 255598bc09dSHong Zhang api[0] = 0; 256598bc09dSHong Zhang 257598bc09dSHong Zhang /* create and initialize a linked list */ 258f84c536bSHong Zhang armax = ad->rmax+ao->rmax; 2599467f45fSHong Zhang if (size >1) { 260f84c536bSHong Zhang prmax = PetscMax(p_loc->rmax,p_oth->rmax); 2619467f45fSHong Zhang } else { 2629467f45fSHong Zhang prmax = p_loc->rmax; 2639467f45fSHong Zhang } 264f84c536bSHong Zhang nlnk_max = armax*prmax; 265f84c536bSHong Zhang if (!nlnk_max || nlnk_max > pN) nlnk_max = pN; 2660d3134e9SHong Zhang ierr = PetscLLCondensedCreate(nlnk_max,pN,&lnk,&lnkbt);CHKERRQ(ierr); 267598bc09dSHong Zhang 268bfcd1a73SHong Zhang /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */ 269bfcd1a73SHong Zhang ierr = PetscFreeSpaceGet((PetscInt)(fill*(adi[am]+aoi[am]+pi_loc[pm])),&free_space);CHKERRQ(ierr); 2702205254eSKarl Rupp 271598bc09dSHong Zhang current_space = free_space; 272598bc09dSHong Zhang 273598bc09dSHong Zhang ierr = MatPreallocateInitialize(comm,am,pn,dnz,onz);CHKERRQ(ierr); 274598bc09dSHong Zhang for (i=0; i<am; i++) { 275598bc09dSHong Zhang /* diagonal portion of A */ 276598bc09dSHong Zhang nzi = adi[i+1] - adi[i]; 277598bc09dSHong Zhang for (j=0; j<nzi; j++) { 278598bc09dSHong Zhang row = *adj++; 279598bc09dSHong Zhang pnz = pi_loc[row+1] - pi_loc[row]; 280598bc09dSHong Zhang Jptr = pj_loc + pi_loc[row]; 281598bc09dSHong Zhang /* add non-zero cols of P into the sorted linked list lnk */ 2821fbbb359SHong Zhang ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 283598bc09dSHong Zhang } 284598bc09dSHong Zhang /* off-diagonal portion of A */ 285598bc09dSHong Zhang nzi = aoi[i+1] - aoi[i]; 286598bc09dSHong Zhang for (j=0; j<nzi; j++) { 287598bc09dSHong Zhang row = *aoj++; 288598bc09dSHong Zhang pnz = pi_oth[row+1] - pi_oth[row]; 289598bc09dSHong Zhang Jptr = pj_oth + pi_oth[row]; 2901fbbb359SHong Zhang ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 291598bc09dSHong Zhang } 292598bc09dSHong Zhang 293d14fa243SHong Zhang apnz = lnk[0]; 294598bc09dSHong Zhang api[i+1] = api[i] + apnz; 295598bc09dSHong Zhang 296598bc09dSHong Zhang /* if free space is not available, double the total space in the list */ 297598bc09dSHong Zhang if (current_space->local_remaining<apnz) { 298598bc09dSHong Zhang ierr = PetscFreeSpaceGet(apnz+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 299598bc09dSHong Zhang nspacedouble++; 300598bc09dSHong Zhang } 301598bc09dSHong Zhang 302598bc09dSHong Zhang /* Copy data into free space, then initialize lnk */ 3031fbbb359SHong Zhang ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 304598bc09dSHong Zhang ierr = MatPreallocateSet(i+rstart,apnz,current_space->array,dnz,onz);CHKERRQ(ierr); 3052205254eSKarl Rupp 306598bc09dSHong Zhang current_space->array += apnz; 307598bc09dSHong Zhang current_space->local_used += apnz; 308598bc09dSHong Zhang current_space->local_remaining -= apnz; 309598bc09dSHong Zhang } 310598bc09dSHong Zhang 311598bc09dSHong Zhang /* Allocate space for apj, initialize apj, and */ 312598bc09dSHong Zhang /* destroy list of free space and other temporary array(s) */ 313854ce69bSBarry Smith ierr = PetscMalloc1(api[am]+1,&ptap->apj);CHKERRQ(ierr); 314a1a4e84aSHong Zhang apj = ptap->apj; 315a1a4e84aSHong Zhang ierr = PetscFreeSpaceContiguous(&free_space,ptap->apj);CHKERRQ(ierr); 316598bc09dSHong Zhang ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 317598bc09dSHong Zhang 318f84c536bSHong Zhang /* malloc apa to store dense row A[i,:]*P */ 3191795a4d1SJed Brown ierr = PetscCalloc1(pN,&apa);CHKERRQ(ierr); 3202205254eSKarl Rupp 321f84c536bSHong Zhang ptap->apa = apa; 322f84c536bSHong Zhang 323bfcd1a73SHong Zhang /* create and assemble symbolic parallel matrix Cmpi */ 324598bc09dSHong Zhang /*----------------------------------------------------*/ 325bfcd1a73SHong Zhang ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr); 326bfcd1a73SHong Zhang ierr = MatSetSizes(Cmpi,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 32733d57670SJed Brown ierr = MatSetBlockSizesFromMats(Cmpi,A,P);CHKERRQ(ierr); 328a2f3521dSMark F. Adams 329bfcd1a73SHong Zhang ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr); 330bfcd1a73SHong Zhang ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr); 331598bc09dSHong Zhang ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 332598bc09dSHong Zhang for (i=0; i<am; i++) { 333598bc09dSHong Zhang row = i + rstart; 334598bc09dSHong Zhang apnz = api[i+1] - api[i]; 335bfcd1a73SHong Zhang ierr = MatSetValues(Cmpi,1,&row,apnz,apj,apa,INSERT_VALUES);CHKERRQ(ierr); 336598bc09dSHong Zhang apj += apnz; 337598bc09dSHong Zhang } 338bfcd1a73SHong Zhang ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 339bfcd1a73SHong Zhang ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 340598bc09dSHong Zhang 341bfcd1a73SHong Zhang ptap->destroy = Cmpi->ops->destroy; 342bfcd1a73SHong Zhang ptap->duplicate = Cmpi->ops->duplicate; 343f8efd71cSHong Zhang Cmpi->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable; 344bfcd1a73SHong Zhang Cmpi->ops->destroy = MatDestroy_MPIAIJ_MatMatMult; 345bfcd1a73SHong Zhang Cmpi->ops->duplicate = MatDuplicate_MPIAIJ_MatMatMult; 346598bc09dSHong Zhang 347bfcd1a73SHong Zhang /* attach the supporting struct to Cmpi for reuse */ 348bfcd1a73SHong Zhang c = (Mat_MPIAIJ*)Cmpi->data; 349598bc09dSHong Zhang c->ptap = ptap; 350598bc09dSHong Zhang 351bfcd1a73SHong Zhang *C = Cmpi; 352bfcd1a73SHong Zhang 353bfcd1a73SHong Zhang /* set MatInfo */ 354118727c9SMark F. Adams afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1) + 1.e-5; 355bfcd1a73SHong Zhang if (afill < 1.0) afill = 1.0; 356bfcd1a73SHong Zhang Cmpi->info.mallocs = nspacedouble; 357bfcd1a73SHong Zhang Cmpi->info.fill_ratio_given = fill; 358bfcd1a73SHong Zhang Cmpi->info.fill_ratio_needed = afill; 359bfcd1a73SHong Zhang 360bfcd1a73SHong Zhang #if defined(PETSC_USE_INFO) 361bfcd1a73SHong Zhang if (api[am]) { 36257622a8eSBarry Smith ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr); 36357622a8eSBarry Smith ierr = PetscInfo1(Cmpi,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr); 364bfcd1a73SHong Zhang } else { 365bfcd1a73SHong Zhang ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr); 366bfcd1a73SHong Zhang } 367bfcd1a73SHong Zhang #endif 368598bc09dSHong Zhang PetscFunctionReturn(0); 369598bc09dSHong Zhang } 370598bc09dSHong Zhang 3719123193aSHong Zhang #undef __FUNCT__ 3729123193aSHong Zhang #define __FUNCT__ "MatMatMult_MPIAIJ_MPIDense" 3739123193aSHong Zhang PetscErrorCode MatMatMult_MPIAIJ_MPIDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 3749123193aSHong Zhang { 3759123193aSHong Zhang PetscErrorCode ierr; 3769123193aSHong Zhang 3779123193aSHong Zhang PetscFunctionBegin; 3789123193aSHong Zhang if (scall == MAT_INITIAL_MATRIX) { 3793ff4c91cSHong Zhang ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr); 3809123193aSHong Zhang ierr = MatMatMultSymbolic_MPIAIJ_MPIDense(A,B,fill,C);CHKERRQ(ierr); 3813ff4c91cSHong Zhang ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr); 3829123193aSHong Zhang } 3833ff4c91cSHong Zhang ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr); 3849123193aSHong Zhang ierr = MatMatMultNumeric_MPIAIJ_MPIDense(A,B,*C);CHKERRQ(ierr); 3853ff4c91cSHong Zhang ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr); 3869123193aSHong Zhang PetscFunctionReturn(0); 3879123193aSHong Zhang } 3889123193aSHong Zhang 3894324174eSBarry Smith typedef struct { 3904324174eSBarry Smith Mat workB; 3914324174eSBarry Smith PetscScalar *rvalues,*svalues; 3924324174eSBarry Smith MPI_Request *rwaits,*swaits; 3934324174eSBarry Smith } MPIAIJ_MPIDense; 3944324174eSBarry Smith 3957af9e4fdSHong Zhang #undef __FUNCT__ 39696b95a6bSBarry Smith #define __FUNCT__ "MatMPIAIJ_MPIDenseDestroy" 39796b95a6bSBarry Smith PetscErrorCode MatMPIAIJ_MPIDenseDestroy(void *ctx) 3984324174eSBarry Smith { 3994324174eSBarry Smith MPIAIJ_MPIDense *contents = (MPIAIJ_MPIDense*) ctx; 4004324174eSBarry Smith PetscErrorCode ierr; 4014324174eSBarry Smith 4024324174eSBarry Smith PetscFunctionBegin; 4036bf464f9SBarry Smith ierr = MatDestroy(&contents->workB);CHKERRQ(ierr); 4044324174eSBarry Smith ierr = PetscFree4(contents->rvalues,contents->svalues,contents->rwaits,contents->swaits);CHKERRQ(ierr); 4054324174eSBarry Smith ierr = PetscFree(contents);CHKERRQ(ierr); 4064324174eSBarry Smith PetscFunctionReturn(0); 4074324174eSBarry Smith } 4084324174eSBarry Smith 4099123193aSHong Zhang #undef __FUNCT__ 410*fd4e9aacSBarry Smith #define __FUNCT__ "MatMatMultNumeric_MPIDense" 411*fd4e9aacSBarry Smith /* 412*fd4e9aacSBarry Smith This is a "dummy function" that handles the case where matrix C was created as a dense matrix 413*fd4e9aacSBarry Smith directly by the user and passed to MatMatMult() with the MAT_REUSE_MATRIX option 414*fd4e9aacSBarry Smith 415*fd4e9aacSBarry Smith It is the same as MatMatMultSymbolic_MPIAIJ_MPIDense() except does not create C 416*fd4e9aacSBarry Smith */ 417*fd4e9aacSBarry Smith PetscErrorCode MatMatMultNumeric_MPIDense(Mat A,Mat B,Mat C) 418*fd4e9aacSBarry Smith { 419*fd4e9aacSBarry Smith PetscErrorCode ierr; 420*fd4e9aacSBarry Smith PetscBool flg; 421*fd4e9aacSBarry Smith Mat_MPIAIJ *aij = (Mat_MPIAIJ*) A->data; 422*fd4e9aacSBarry Smith PetscInt nz = aij->B->cmap->n; 423*fd4e9aacSBarry Smith PetscContainer container; 424*fd4e9aacSBarry Smith MPIAIJ_MPIDense *contents; 425*fd4e9aacSBarry Smith VecScatter ctx = aij->Mvctx; 426*fd4e9aacSBarry Smith VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata; 427*fd4e9aacSBarry Smith VecScatter_MPI_General *to = (VecScatter_MPI_General*) ctx->todata; 428*fd4e9aacSBarry Smith 429*fd4e9aacSBarry Smith PetscFunctionBegin; 430*fd4e9aacSBarry Smith ierr = PetscObjectTypeCompare((PetscObject)B,MATMPIDENSE,&flg);CHKERRQ(ierr); 431*fd4e9aacSBarry Smith if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Second matrix must be mpidense"); 432*fd4e9aacSBarry Smith 433*fd4e9aacSBarry Smith /* Handle case where where user provided the final C matrix rather than calling MatMatMult() with MAT_INITIAL_MATRIX*/ 434*fd4e9aacSBarry Smith ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&flg);CHKERRQ(ierr); 435*fd4e9aacSBarry Smith if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"First matrix must be MPIAIJ"); 436*fd4e9aacSBarry Smith 437*fd4e9aacSBarry Smith C->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIDense; 438*fd4e9aacSBarry Smith 439*fd4e9aacSBarry Smith ierr = PetscNew(&contents);CHKERRQ(ierr); 440*fd4e9aacSBarry Smith /* Create work matrix used to store off processor rows of B needed for local product */ 441*fd4e9aacSBarry Smith ierr = MatCreateSeqDense(PETSC_COMM_SELF,nz,B->cmap->N,NULL,&contents->workB);CHKERRQ(ierr); 442*fd4e9aacSBarry Smith /* Create work arrays needed */ 443*fd4e9aacSBarry Smith ierr = PetscMalloc4(B->cmap->N*from->starts[from->n],&contents->rvalues, 444*fd4e9aacSBarry Smith B->cmap->N*to->starts[to->n],&contents->svalues, 445*fd4e9aacSBarry Smith from->n,&contents->rwaits, 446*fd4e9aacSBarry Smith to->n,&contents->swaits);CHKERRQ(ierr); 447*fd4e9aacSBarry Smith 448*fd4e9aacSBarry Smith ierr = PetscContainerCreate(PetscObjectComm((PetscObject)A),&container);CHKERRQ(ierr); 449*fd4e9aacSBarry Smith ierr = PetscContainerSetPointer(container,contents);CHKERRQ(ierr); 450*fd4e9aacSBarry Smith ierr = PetscContainerSetUserDestroy(container,MatMPIAIJ_MPIDenseDestroy);CHKERRQ(ierr); 451*fd4e9aacSBarry Smith ierr = PetscObjectCompose((PetscObject)C,"workB",(PetscObject)container);CHKERRQ(ierr); 452*fd4e9aacSBarry Smith ierr = PetscContainerDestroy(&container);CHKERRQ(ierr); 453*fd4e9aacSBarry Smith 454*fd4e9aacSBarry Smith ierr = (*C->ops->matmultnumeric)(A,B,C);CHKERRQ(ierr); 455*fd4e9aacSBarry Smith PetscFunctionReturn(0); 456*fd4e9aacSBarry Smith } 457*fd4e9aacSBarry Smith 458*fd4e9aacSBarry Smith #undef __FUNCT__ 4599123193aSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIDense" 4609123193aSHong Zhang PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIDense(Mat A,Mat B,PetscReal fill,Mat *C) 4619123193aSHong Zhang { 4629123193aSHong Zhang PetscErrorCode ierr; 463f32d5d43SBarry Smith Mat_MPIAIJ *aij = (Mat_MPIAIJ*) A->data; 464d0f46423SBarry Smith PetscInt nz = aij->B->cmap->n; 465bf0cc555SLisandro Dalcin PetscContainer container; 4664324174eSBarry Smith MPIAIJ_MPIDense *contents; 4674324174eSBarry Smith VecScatter ctx = aij->Mvctx; 4684324174eSBarry Smith VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata; 4694324174eSBarry Smith VecScatter_MPI_General *to = (VecScatter_MPI_General*) ctx->todata; 470d0f46423SBarry Smith PetscInt m = A->rmap->n,n=B->cmap->n; 4719123193aSHong Zhang 4729123193aSHong Zhang PetscFunctionBegin; 473ce94432eSBarry Smith ierr = MatCreate(PetscObjectComm((PetscObject)B),C);CHKERRQ(ierr); 474d0f46423SBarry Smith ierr = MatSetSizes(*C,m,n,A->rmap->N,B->cmap->N);CHKERRQ(ierr); 47533d57670SJed Brown ierr = MatSetBlockSizesFromMats(*C,A,B);CHKERRQ(ierr); 476cb2480eaSBarry Smith ierr = MatSetType(*C,MATMPIDENSE);CHKERRQ(ierr); 4770298fd71SBarry Smith ierr = MatMPIDenseSetPreallocation(*C,NULL);CHKERRQ(ierr); 478cb2480eaSBarry Smith ierr = MatAssemblyBegin(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 479cb2480eaSBarry Smith ierr = MatAssemblyEnd(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 4802205254eSKarl Rupp 481d73949e8SHong Zhang (*C)->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIDense; 482f32d5d43SBarry Smith 483b00a9115SJed Brown ierr = PetscNew(&contents);CHKERRQ(ierr); 484f32d5d43SBarry Smith /* Create work matrix used to store off processor rows of B needed for local product */ 4850298fd71SBarry Smith ierr = MatCreateSeqDense(PETSC_COMM_SELF,nz,B->cmap->N,NULL,&contents->workB);CHKERRQ(ierr); 4864324174eSBarry Smith /* Create work arrays needed */ 487dcca6d9dSJed Brown ierr = PetscMalloc4(B->cmap->N*from->starts[from->n],&contents->rvalues, 488dcca6d9dSJed Brown B->cmap->N*to->starts[to->n],&contents->svalues, 489dcca6d9dSJed Brown from->n,&contents->rwaits, 490dcca6d9dSJed Brown to->n,&contents->swaits);CHKERRQ(ierr); 4914324174eSBarry Smith 492ce94432eSBarry Smith ierr = PetscContainerCreate(PetscObjectComm((PetscObject)A),&container);CHKERRQ(ierr); 493bf0cc555SLisandro Dalcin ierr = PetscContainerSetPointer(container,contents);CHKERRQ(ierr); 49496b95a6bSBarry Smith ierr = PetscContainerSetUserDestroy(container,MatMPIAIJ_MPIDenseDestroy);CHKERRQ(ierr); 495bf0cc555SLisandro Dalcin ierr = PetscObjectCompose((PetscObject)(*C),"workB",(PetscObject)container);CHKERRQ(ierr); 496bf0cc555SLisandro Dalcin ierr = PetscContainerDestroy(&container);CHKERRQ(ierr); 4979123193aSHong Zhang PetscFunctionReturn(0); 4989123193aSHong Zhang } 4999123193aSHong Zhang 5007af9e4fdSHong Zhang #undef __FUNCT__ 5017af9e4fdSHong Zhang #define __FUNCT__ "MatMPIDenseScatter" 502f32d5d43SBarry Smith /* 5032636ff26SBarry Smith Performs an efficient scatter on the rows of B needed by this process; this is 5042636ff26SBarry Smith a modification of the VecScatterBegin_() routines. 505f32d5d43SBarry Smith */ 5064324174eSBarry Smith PetscErrorCode MatMPIDenseScatter(Mat A,Mat B,Mat C,Mat *outworkB) 507f32d5d43SBarry Smith { 508f32d5d43SBarry Smith Mat_MPIAIJ *aij = (Mat_MPIAIJ*)A->data; 509f32d5d43SBarry Smith PetscErrorCode ierr; 510f32d5d43SBarry Smith PetscScalar *b,*w,*svalues,*rvalues; 511f32d5d43SBarry Smith VecScatter ctx = aij->Mvctx; 512f32d5d43SBarry Smith VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata; 513f32d5d43SBarry Smith VecScatter_MPI_General *to = (VecScatter_MPI_General*) ctx->todata; 514f32d5d43SBarry Smith PetscInt i,j,k; 515aa5bb8c0SSatish Balay PetscInt *sindices,*sstarts,*rindices,*rstarts; 516aa5bb8c0SSatish Balay PetscMPIInt *sprocs,*rprocs,nrecvs; 517f32d5d43SBarry Smith MPI_Request *swaits,*rwaits; 518ce94432eSBarry Smith MPI_Comm comm; 519d0f46423SBarry Smith PetscMPIInt tag = ((PetscObject)ctx)->tag,ncols = B->cmap->N, nrows = aij->B->cmap->n,imdex,nrowsB = B->rmap->n; 520f32d5d43SBarry Smith MPI_Status status; 5214324174eSBarry Smith MPIAIJ_MPIDense *contents; 522bf0cc555SLisandro Dalcin PetscContainer container; 5234324174eSBarry Smith Mat workB; 524f32d5d43SBarry Smith 525f32d5d43SBarry Smith PetscFunctionBegin; 526ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 527bf0cc555SLisandro Dalcin ierr = PetscObjectQuery((PetscObject)C,"workB",(PetscObject*)&container);CHKERRQ(ierr); 52829825010SHong Zhang if (!container) SETERRQ(comm,PETSC_ERR_PLIB,"Container does not exist"); 529bf0cc555SLisandro Dalcin ierr = PetscContainerGetPointer(container,(void**)&contents);CHKERRQ(ierr); 5304324174eSBarry Smith 5314324174eSBarry Smith workB = *outworkB = contents->workB; 532cf1a0672SHong 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); 533f32d5d43SBarry Smith sindices = to->indices; 534f32d5d43SBarry Smith sstarts = to->starts; 535f32d5d43SBarry Smith sprocs = to->procs; 5364324174eSBarry Smith swaits = contents->swaits; 5374324174eSBarry Smith svalues = contents->svalues; 538f32d5d43SBarry Smith 539f32d5d43SBarry Smith rindices = from->indices; 540f32d5d43SBarry Smith rstarts = from->starts; 541f32d5d43SBarry Smith rprocs = from->procs; 5424324174eSBarry Smith rwaits = contents->rwaits; 5434324174eSBarry Smith rvalues = contents->rvalues; 544f32d5d43SBarry Smith 5458c778c55SBarry Smith ierr = MatDenseGetArray(B,&b);CHKERRQ(ierr); 5468c778c55SBarry Smith ierr = MatDenseGetArray(workB,&w);CHKERRQ(ierr); 547f32d5d43SBarry Smith 548f32d5d43SBarry Smith for (i=0; i<from->n; i++) { 549f32d5d43SBarry Smith ierr = MPI_Irecv(rvalues+ncols*rstarts[i],ncols*(rstarts[i+1]-rstarts[i]),MPIU_SCALAR,rprocs[i],tag,comm,rwaits+i);CHKERRQ(ierr); 550f32d5d43SBarry Smith } 551f32d5d43SBarry Smith 552f32d5d43SBarry Smith for (i=0; i<to->n; i++) { 553f32d5d43SBarry Smith /* pack a message at a time */ 554f32d5d43SBarry Smith for (j=0; j<sstarts[i+1]-sstarts[i]; j++) { 555f32d5d43SBarry Smith for (k=0; k<ncols; k++) { 556f32d5d43SBarry Smith svalues[ncols*(sstarts[i] + j) + k] = b[sindices[sstarts[i]+j] + nrowsB*k]; 5572636ff26SBarry Smith } 5582636ff26SBarry Smith } 559f32d5d43SBarry Smith ierr = MPI_Isend(svalues+ncols*sstarts[i],ncols*(sstarts[i+1]-sstarts[i]),MPIU_SCALAR,sprocs[i],tag,comm,swaits+i);CHKERRQ(ierr); 560f32d5d43SBarry Smith } 5612636ff26SBarry Smith 562f32d5d43SBarry Smith nrecvs = from->n; 563f32d5d43SBarry Smith while (nrecvs) { 564f32d5d43SBarry Smith ierr = MPI_Waitany(from->n,rwaits,&imdex,&status);CHKERRQ(ierr); 565f32d5d43SBarry Smith nrecvs--; 566f32d5d43SBarry Smith /* unpack a message at a time */ 567f32d5d43SBarry Smith for (j=0; j<rstarts[imdex+1]-rstarts[imdex]; j++) { 568f32d5d43SBarry Smith for (k=0; k<ncols; k++) { 569f32d5d43SBarry Smith w[rindices[rstarts[imdex]+j] + nrows*k] = rvalues[ncols*(rstarts[imdex] + j) + k]; 5702636ff26SBarry Smith } 5712636ff26SBarry Smith } 572f32d5d43SBarry Smith } 573cb9801acSJed Brown if (to->n) {ierr = MPI_Waitall(to->n,swaits,to->sstatus);CHKERRQ(ierr);} 574f32d5d43SBarry Smith 5758c778c55SBarry Smith ierr = MatDenseRestoreArray(B,&b);CHKERRQ(ierr); 5768c778c55SBarry Smith ierr = MatDenseRestoreArray(workB,&w);CHKERRQ(ierr); 577f32d5d43SBarry Smith ierr = MatAssemblyBegin(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 578f32d5d43SBarry Smith ierr = MatAssemblyEnd(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 579f32d5d43SBarry Smith PetscFunctionReturn(0); 580f32d5d43SBarry Smith } 581f32d5d43SBarry Smith extern PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat,Mat,Mat); 582f32d5d43SBarry Smith 5839123193aSHong Zhang #undef __FUNCT__ 5849123193aSHong Zhang #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIDense" 5859123193aSHong Zhang PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIDense(Mat A,Mat B,Mat C) 5869123193aSHong Zhang { 5879123193aSHong Zhang PetscErrorCode ierr; 588f32d5d43SBarry Smith Mat_MPIAIJ *aij = (Mat_MPIAIJ*)A->data; 589f32d5d43SBarry Smith Mat_MPIDense *bdense = (Mat_MPIDense*)B->data; 590f32d5d43SBarry Smith Mat_MPIDense *cdense = (Mat_MPIDense*)C->data; 591f32d5d43SBarry Smith Mat workB; 5929123193aSHong Zhang 5939123193aSHong Zhang PetscFunctionBegin; 594f32d5d43SBarry Smith /* diagonal block of A times all local rows of B*/ 595f32d5d43SBarry Smith ierr = MatMatMultNumeric_SeqAIJ_SeqDense(aij->A,bdense->A,cdense->A);CHKERRQ(ierr); 596f32d5d43SBarry Smith 597f32d5d43SBarry Smith /* get off processor parts of B needed to complete the product */ 5984324174eSBarry Smith ierr = MatMPIDenseScatter(A,B,C,&workB);CHKERRQ(ierr); 599f32d5d43SBarry Smith 600f32d5d43SBarry Smith /* off-diagonal block of A times nonlocal rows of B */ 601f32d5d43SBarry Smith ierr = MatMatMultNumericAdd_SeqAIJ_SeqDense(aij->B,workB,cdense->A);CHKERRQ(ierr); 6029123193aSHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 6039123193aSHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 6049123193aSHong Zhang PetscFunctionReturn(0); 6059123193aSHong Zhang } 606cf1a0672SHong Zhang 6071634b032SHong Zhang #undef __FUNCT__ 608f8efd71cSHong Zhang #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIAIJ" 609f8efd71cSHong Zhang PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C) 6101634b032SHong Zhang { 6111634b032SHong Zhang PetscErrorCode ierr; 612cf1a0672SHong Zhang Mat_MPIAIJ *a = (Mat_MPIAIJ*)A->data,*c=(Mat_MPIAIJ*)C->data; 613cf1a0672SHong Zhang Mat_SeqAIJ *ad = (Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data; 614cf1a0672SHong Zhang Mat_SeqAIJ *cd = (Mat_SeqAIJ*)(c->A)->data,*co=(Mat_SeqAIJ*)(c->B)->data; 615cf1a0672SHong Zhang PetscInt *adi = ad->i,*adj,*aoi=ao->i,*aoj; 616cf1a0672SHong Zhang PetscScalar *ada,*aoa,*cda=cd->a,*coa=co->a; 617cf1a0672SHong Zhang Mat_SeqAIJ *p_loc,*p_oth; 618cf1a0672SHong Zhang PetscInt *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pj; 61929825010SHong Zhang PetscScalar *pa_loc,*pa_oth,*pa,valtmp,*ca; 62029825010SHong Zhang PetscInt cm = C->rmap->n,anz,pnz; 621cf1a0672SHong Zhang Mat_PtAPMPI *ptap = c->ptap; 62229825010SHong Zhang PetscScalar *apa_sparse = ptap->apa; 623cf1a0672SHong Zhang PetscInt *api,*apj,*apJ,i,j,k,row; 62429825010SHong Zhang PetscInt cstart = C->cmap->rstart; 62529825010SHong Zhang PetscInt cdnz,conz,k0,k1,nextp; 626a7c7454dSHong Zhang MPI_Comm comm; 627a7c7454dSHong Zhang PetscMPIInt size; 6281634b032SHong Zhang 6291634b032SHong Zhang PetscFunctionBegin; 630a7c7454dSHong Zhang ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 631a7c7454dSHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 632a7c7454dSHong Zhang 633cf1a0672SHong Zhang /* 1) get P_oth = ptap->P_oth and P_loc = ptap->P_loc */ 634cf1a0672SHong Zhang /*-----------------------------------------------------*/ 635cf1a0672SHong Zhang /* update numerical values of P_oth and P_loc */ 636b7f45c76SHong Zhang ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 637cf1a0672SHong Zhang ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 638cf1a0672SHong Zhang 639cf1a0672SHong Zhang /* 2) compute numeric C_loc = A_loc*P = Ad*P_loc + Ao*P_oth */ 640cf1a0672SHong Zhang /*----------------------------------------------------------*/ 641cf1a0672SHong Zhang /* get data from symbolic products */ 642cf1a0672SHong Zhang p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 643cf1a0672SHong Zhang pi_loc = p_loc->i; pj_loc = p_loc->j; pa_loc = p_loc->a; 644a7c7454dSHong Zhang if (size >1) { 645a7c7454dSHong Zhang p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 646cf1a0672SHong Zhang pi_oth = p_oth->i; pj_oth = p_oth->j; pa_oth = p_oth->a; 64720e3dc75SHong Zhang } else { 64820e3dc75SHong Zhang p_oth = NULL; pi_oth = NULL; pj_oth = NULL; pa_oth = NULL; 649a7c7454dSHong Zhang } 650cf1a0672SHong Zhang 65129825010SHong Zhang api = ptap->api; 65229825010SHong Zhang apj = ptap->apj; 653cf1a0672SHong Zhang for (i=0; i<cm; i++) { 65429825010SHong Zhang apJ = apj + api[i]; 65529825010SHong Zhang 656cf1a0672SHong Zhang /* diagonal portion of A */ 657cf1a0672SHong Zhang anz = adi[i+1] - adi[i]; 658cf1a0672SHong Zhang adj = ad->j + adi[i]; 659cf1a0672SHong Zhang ada = ad->a + adi[i]; 660cf1a0672SHong Zhang for (j=0; j<anz; j++) { 661cf1a0672SHong Zhang row = adj[j]; 662cf1a0672SHong Zhang pnz = pi_loc[row+1] - pi_loc[row]; 663cf1a0672SHong Zhang pj = pj_loc + pi_loc[row]; 664cf1a0672SHong Zhang pa = pa_loc + pi_loc[row]; 66529825010SHong Zhang /* perform sparse axpy */ 666cf1a0672SHong Zhang valtmp = ada[j]; 66729825010SHong Zhang nextp = 0; 66829825010SHong Zhang for (k=0; nextp<pnz; k++) { 66929825010SHong Zhang if (apJ[k] == pj[nextp]) { /* column of AP == column of P */ 67029825010SHong Zhang apa_sparse[k] += valtmp*pa[nextp++]; 67129825010SHong Zhang } 6721634b032SHong Zhang } 673cf1a0672SHong Zhang ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr); 674cf1a0672SHong Zhang } 6751634b032SHong Zhang 676cf1a0672SHong Zhang /* off-diagonal portion of A */ 677cf1a0672SHong Zhang anz = aoi[i+1] - aoi[i]; 678cf1a0672SHong Zhang aoj = ao->j + aoi[i]; 679cf1a0672SHong Zhang aoa = ao->a + aoi[i]; 680cf1a0672SHong Zhang for (j=0; j<anz; j++) { 681cf1a0672SHong Zhang row = aoj[j]; 682cf1a0672SHong Zhang pnz = pi_oth[row+1] - pi_oth[row]; 683cf1a0672SHong Zhang pj = pj_oth + pi_oth[row]; 684cf1a0672SHong Zhang pa = pa_oth + pi_oth[row]; 68529825010SHong Zhang /* perform sparse axpy */ 686cf1a0672SHong Zhang valtmp = aoa[j]; 68729825010SHong Zhang nextp = 0; 68829825010SHong Zhang for (k=0; nextp<pnz; k++) { 68929825010SHong Zhang if (apJ[k] == pj[nextp]) { /* column of AP == column of P */ 69029825010SHong Zhang apa_sparse[k] += valtmp*pa[nextp++]; 69129825010SHong Zhang } 692cf1a0672SHong Zhang } 693cf1a0672SHong Zhang ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr); 694cf1a0672SHong Zhang } 6951634b032SHong Zhang 696cf1a0672SHong Zhang /* set values in C */ 697cf1a0672SHong Zhang cdnz = cd->i[i+1] - cd->i[i]; 698cf1a0672SHong Zhang conz = co->i[i+1] - co->i[i]; 6991634b032SHong Zhang 700cf1a0672SHong Zhang /* 1st off-diagoanl part of C */ 701cf1a0672SHong Zhang ca = coa + co->i[i]; 702cf1a0672SHong Zhang k = 0; 703cf1a0672SHong Zhang for (k0=0; k0<conz; k0++) { 704cf1a0672SHong Zhang if (apJ[k] >= cstart) break; 70529825010SHong Zhang ca[k0] = apa_sparse[k]; 70629825010SHong Zhang apa_sparse[k] = 0.0; 707cf1a0672SHong Zhang k++; 708cf1a0672SHong Zhang } 7091634b032SHong Zhang 710cf1a0672SHong Zhang /* diagonal part of C */ 711cf1a0672SHong Zhang ca = cda + cd->i[i]; 712cf1a0672SHong Zhang for (k1=0; k1<cdnz; k1++) { 71329825010SHong Zhang ca[k1] = apa_sparse[k]; 71429825010SHong Zhang apa_sparse[k] = 0.0; 715cf1a0672SHong Zhang k++; 716cf1a0672SHong Zhang } 717cf1a0672SHong Zhang 718cf1a0672SHong Zhang /* 2nd off-diagoanl part of C */ 719cf1a0672SHong Zhang ca = coa + co->i[i]; 720cf1a0672SHong Zhang for (; k0<conz; k0++) { 72129825010SHong Zhang ca[k0] = apa_sparse[k]; 72229825010SHong Zhang apa_sparse[k] = 0.0; 723cf1a0672SHong Zhang k++; 724cf1a0672SHong Zhang } 725cf1a0672SHong Zhang } 726cf1a0672SHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 727cf1a0672SHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 728cf1a0672SHong Zhang PetscFunctionReturn(0); 729cf1a0672SHong Zhang } 730cf1a0672SHong Zhang 7310fc8cf34SHong Zhang /* same as MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(), except using LLCondensed to avoid O(BN) memory requirement */ 732cf1a0672SHong Zhang #undef __FUNCT__ 733b2405163SHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ" 734b2405163SHong Zhang PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C) 735cf1a0672SHong Zhang { 736cf1a0672SHong Zhang PetscErrorCode ierr; 737ce94432eSBarry Smith MPI_Comm comm; 738a7c7454dSHong Zhang PetscMPIInt size; 739cf1a0672SHong Zhang Mat Cmpi; 740cf1a0672SHong Zhang Mat_PtAPMPI *ptap; 7410298fd71SBarry Smith PetscFreeSpaceList free_space = NULL,current_space=NULL; 742cf1a0672SHong Zhang Mat_MPIAIJ *a = (Mat_MPIAIJ*)A->data,*c; 743cf1a0672SHong Zhang Mat_SeqAIJ *ad = (Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_loc,*p_oth; 744cf1a0672SHong Zhang PetscInt *pi_loc,*pj_loc,*pi_oth,*pj_oth,*dnz,*onz; 745cf1a0672SHong Zhang PetscInt *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,rstart=A->rmap->rstart; 746f84c536bSHong Zhang PetscInt i,pnz,row,*api,*apj,*Jptr,apnz,nspacedouble=0,j,nzi,*lnk,apnz_max=0; 747cf1a0672SHong Zhang PetscInt am=A->rmap->n,pN=P->cmap->N,pn=P->cmap->n,pm=P->rmap->n; 748cf1a0672SHong Zhang PetscInt nlnk_max,armax,prmax; 749cf1a0672SHong Zhang PetscReal afill; 750cf1a0672SHong Zhang PetscScalar *apa; 751cf1a0672SHong Zhang 752cf1a0672SHong Zhang PetscFunctionBegin; 753ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 754a7c7454dSHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 755a7c7454dSHong Zhang 756cf1a0672SHong Zhang /* create struct Mat_PtAPMPI and attached it to C later */ 757b00a9115SJed Brown ierr = PetscNew(&ptap);CHKERRQ(ierr); 758cf1a0672SHong Zhang 759cf1a0672SHong Zhang /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */ 760b7f45c76SHong Zhang ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 76119f0e57aSHong Zhang 762cf1a0672SHong Zhang /* get P_loc by taking all local rows of P */ 763cf1a0672SHong Zhang ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 764cf1a0672SHong Zhang 765cf1a0672SHong Zhang p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 766cf1a0672SHong Zhang pi_loc = p_loc->i; pj_loc = p_loc->j; 767a7c7454dSHong Zhang if (size > 1) { 768a7c7454dSHong Zhang p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 76920e3dc75SHong Zhang pi_oth = p_oth->i; pj_oth = p_oth->j; 770968382fdSHong Zhang } else { 77120e3dc75SHong Zhang p_oth = NULL; 77220e3dc75SHong Zhang pi_oth = NULL; pj_oth = NULL; 773a7c7454dSHong Zhang } 774cf1a0672SHong Zhang 775cf1a0672SHong Zhang /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */ 776cf1a0672SHong Zhang /*-------------------------------------------------------------------*/ 777854ce69bSBarry Smith ierr = PetscMalloc1(am+2,&api);CHKERRQ(ierr); 778cf1a0672SHong Zhang ptap->api = api; 779cf1a0672SHong Zhang api[0] = 0; 780cf1a0672SHong Zhang 781cf1a0672SHong Zhang /* create and initialize a linked list */ 782cf1a0672SHong Zhang armax = ad->rmax+ao->rmax; 783a7c7454dSHong Zhang if (size >1) { 784cf1a0672SHong Zhang prmax = PetscMax(p_loc->rmax,p_oth->rmax); 785a7c7454dSHong Zhang } else { 786a7c7454dSHong Zhang prmax = p_loc->rmax; 787a7c7454dSHong Zhang } 788cf1a0672SHong Zhang nlnk_max = armax*prmax; 789cf1a0672SHong Zhang if (!nlnk_max || nlnk_max > pN) nlnk_max = pN; 790f84c536bSHong Zhang ierr = PetscLLCondensedCreate_Scalable(nlnk_max,&lnk);CHKERRQ(ierr); 791cf1a0672SHong Zhang 792cf1a0672SHong Zhang /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */ 793cf1a0672SHong Zhang ierr = PetscFreeSpaceGet((PetscInt)(fill*(adi[am]+aoi[am]+pi_loc[pm])),&free_space);CHKERRQ(ierr); 7942205254eSKarl Rupp 795cf1a0672SHong Zhang current_space = free_space; 796cf1a0672SHong Zhang 797cf1a0672SHong Zhang ierr = MatPreallocateInitialize(comm,am,pn,dnz,onz);CHKERRQ(ierr); 798cf1a0672SHong Zhang for (i=0; i<am; i++) { 799cf1a0672SHong Zhang /* diagonal portion of A */ 800cf1a0672SHong Zhang nzi = adi[i+1] - adi[i]; 801cf1a0672SHong Zhang for (j=0; j<nzi; j++) { 802cf1a0672SHong Zhang row = *adj++; 803cf1a0672SHong Zhang pnz = pi_loc[row+1] - pi_loc[row]; 804cf1a0672SHong Zhang Jptr = pj_loc + pi_loc[row]; 805cf1a0672SHong Zhang /* add non-zero cols of P into the sorted linked list lnk */ 806f84c536bSHong Zhang ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr); 807cf1a0672SHong Zhang } 808cf1a0672SHong Zhang /* off-diagonal portion of A */ 809cf1a0672SHong Zhang nzi = aoi[i+1] - aoi[i]; 810cf1a0672SHong Zhang for (j=0; j<nzi; j++) { 811cf1a0672SHong Zhang row = *aoj++; 812cf1a0672SHong Zhang pnz = pi_oth[row+1] - pi_oth[row]; 813cf1a0672SHong Zhang Jptr = pj_oth + pi_oth[row]; 814f84c536bSHong Zhang ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr); 815cf1a0672SHong Zhang } 816cf1a0672SHong Zhang 817f84c536bSHong Zhang apnz = *lnk; 818cf1a0672SHong Zhang api[i+1] = api[i] + apnz; 819cf1a0672SHong Zhang if (apnz > apnz_max) apnz_max = apnz; 820cf1a0672SHong Zhang 821cf1a0672SHong Zhang /* if free space is not available, double the total space in the list */ 822cf1a0672SHong Zhang if (current_space->local_remaining<apnz) { 823cf1a0672SHong Zhang ierr = PetscFreeSpaceGet(apnz+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 824cf1a0672SHong Zhang nspacedouble++; 825cf1a0672SHong Zhang } 826cf1a0672SHong Zhang 827cf1a0672SHong Zhang /* Copy data into free space, then initialize lnk */ 828f84c536bSHong Zhang ierr = PetscLLCondensedClean_Scalable(apnz,current_space->array,lnk);CHKERRQ(ierr); 829cf1a0672SHong Zhang ierr = MatPreallocateSet(i+rstart,apnz,current_space->array,dnz,onz);CHKERRQ(ierr); 8302205254eSKarl Rupp 831cf1a0672SHong Zhang current_space->array += apnz; 832cf1a0672SHong Zhang current_space->local_used += apnz; 833cf1a0672SHong Zhang current_space->local_remaining -= apnz; 834cf1a0672SHong Zhang } 835cf1a0672SHong Zhang 836cf1a0672SHong Zhang /* Allocate space for apj, initialize apj, and */ 837cf1a0672SHong Zhang /* destroy list of free space and other temporary array(s) */ 838854ce69bSBarry Smith ierr = PetscMalloc1(api[am]+1,&ptap->apj);CHKERRQ(ierr); 839cf1a0672SHong Zhang apj = ptap->apj; 840cf1a0672SHong Zhang ierr = PetscFreeSpaceContiguous(&free_space,ptap->apj);CHKERRQ(ierr); 841f84c536bSHong Zhang ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr); 842cf1a0672SHong Zhang 843cf1a0672SHong Zhang /* create and assemble symbolic parallel matrix Cmpi */ 844cf1a0672SHong Zhang /*----------------------------------------------------*/ 845cf1a0672SHong Zhang ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr); 846cf1a0672SHong Zhang ierr = MatSetSizes(Cmpi,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 84733d57670SJed Brown ierr = MatSetBlockSizesFromMats(Cmpi,A,P);CHKERRQ(ierr); 848cf1a0672SHong Zhang ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr); 849cf1a0672SHong Zhang ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr); 850cf1a0672SHong Zhang ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 851cf1a0672SHong Zhang 85229825010SHong Zhang /* malloc apa for assembly Cmpi */ 8531795a4d1SJed Brown ierr = PetscCalloc1(apnz_max,&apa);CHKERRQ(ierr); 8542205254eSKarl Rupp 85529825010SHong Zhang ptap->apa = apa; 856cf1a0672SHong Zhang for (i=0; i<am; i++) { 857cf1a0672SHong Zhang row = i + rstart; 858cf1a0672SHong Zhang apnz = api[i+1] - api[i]; 859cf1a0672SHong Zhang ierr = MatSetValues(Cmpi,1,&row,apnz,apj,apa,INSERT_VALUES);CHKERRQ(ierr); 860cf1a0672SHong Zhang apj += apnz; 861cf1a0672SHong Zhang } 862cf1a0672SHong Zhang ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 863cf1a0672SHong Zhang ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 864cf1a0672SHong Zhang 865cf1a0672SHong Zhang ptap->destroy = Cmpi->ops->destroy; 866cf1a0672SHong Zhang ptap->duplicate = Cmpi->ops->duplicate; 867f8efd71cSHong Zhang Cmpi->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ; 868cf1a0672SHong Zhang Cmpi->ops->destroy = MatDestroy_MPIAIJ_MatMatMult; 869cf1a0672SHong Zhang Cmpi->ops->duplicate = MatDuplicate_MPIAIJ_MatMatMult; 870cf1a0672SHong Zhang 871cf1a0672SHong Zhang /* attach the supporting struct to Cmpi for reuse */ 872cf1a0672SHong Zhang c = (Mat_MPIAIJ*)Cmpi->data; 873cf1a0672SHong Zhang c->ptap = ptap; 874cf1a0672SHong Zhang 875cf1a0672SHong Zhang *C = Cmpi; 876cf1a0672SHong Zhang 877cf1a0672SHong Zhang /* set MatInfo */ 878118727c9SMark F. Adams afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1) + 1.e-5; 879cf1a0672SHong Zhang if (afill < 1.0) afill = 1.0; 880cf1a0672SHong Zhang Cmpi->info.mallocs = nspacedouble; 881cf1a0672SHong Zhang Cmpi->info.fill_ratio_given = fill; 882cf1a0672SHong Zhang Cmpi->info.fill_ratio_needed = afill; 883cf1a0672SHong Zhang 884cf1a0672SHong Zhang #if defined(PETSC_USE_INFO) 885cf1a0672SHong Zhang if (api[am]) { 88657622a8eSBarry Smith ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr); 88757622a8eSBarry Smith ierr = PetscInfo1(Cmpi,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr); 888cf1a0672SHong Zhang } else { 889cf1a0672SHong Zhang ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr); 890cf1a0672SHong Zhang } 891cf1a0672SHong Zhang #endif 8921634b032SHong Zhang PetscFunctionReturn(0); 8931634b032SHong Zhang } 894f96d37fcSHong Zhang 895f96d37fcSHong Zhang /*-------------------------------------------------------------------------*/ 896f96d37fcSHong Zhang #undef __FUNCT__ 897f96d37fcSHong Zhang #define __FUNCT__ "MatTransposeMatMult_MPIAIJ_MPIAIJ" 898c5af039cSHong Zhang PetscErrorCode MatTransposeMatMult_MPIAIJ_MPIAIJ(Mat P,Mat A,MatReuse scall,PetscReal fill,Mat *C) 899f96d37fcSHong Zhang { 900f96d37fcSHong Zhang PetscErrorCode ierr; 901c216dbf3SHong Zhang const char *algTypes[3] = {"scalable","nonscalable","matmatmult"}; 902c216dbf3SHong Zhang PetscInt alg=0; /* set default algorithm */ 903f96d37fcSHong Zhang 904f96d37fcSHong Zhang PetscFunctionBegin; 905c216dbf3SHong Zhang if (scall == MAT_INITIAL_MATRIX) { 906c216dbf3SHong Zhang ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr); 907c216dbf3SHong Zhang ierr = PetscOptionsEList("-mattransposematmult_via","Algorithmic approach","MatTransposeMatMult",algTypes,3,algTypes[0],&alg,NULL);CHKERRQ(ierr); 908c216dbf3SHong Zhang ierr = PetscOptionsEnd();CHKERRQ(ierr); 909c216dbf3SHong Zhang 910c216dbf3SHong Zhang ierr = PetscLogEventBegin(MAT_TransposeMatMultSymbolic,P,A,0,0);CHKERRQ(ierr); 911c216dbf3SHong Zhang switch (alg) { 912c216dbf3SHong Zhang case 1: 9132bbb1c24SHong Zhang ierr = MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(P,A,fill,C);CHKERRQ(ierr); 914c216dbf3SHong Zhang break; 915c216dbf3SHong Zhang case 2: 916275476c6SMatthew G. Knepley { 917ab1b013aSHong Zhang Mat Pt; 918ab1b013aSHong Zhang Mat_PtAPMPI *ptap; 919ab1b013aSHong Zhang Mat_MPIAIJ *c; 920ab1b013aSHong Zhang ierr = MatTranspose(P,MAT_INITIAL_MATRIX,&Pt);CHKERRQ(ierr); 921ab1b013aSHong Zhang ierr = MatMatMult(Pt,A,MAT_INITIAL_MATRIX,fill,C);CHKERRQ(ierr); 922ab1b013aSHong Zhang c = (Mat_MPIAIJ*)(*C)->data; 923ab1b013aSHong Zhang ptap = c->ptap; 924ab1b013aSHong Zhang ptap->Pt = Pt; 925c216dbf3SHong Zhang (*C)->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult; 926c216dbf3SHong Zhang PetscFunctionReturn(0); 927275476c6SMatthew G. Knepley } 928c216dbf3SHong Zhang break; 929c216dbf3SHong Zhang default: 9306da69ca6SHong Zhang ierr = MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(P,A,fill,C);CHKERRQ(ierr); 931c216dbf3SHong Zhang break; 932c216dbf3SHong Zhang } 933c216dbf3SHong Zhang ierr = PetscLogEventEnd(MAT_TransposeMatMultSymbolic,P,A,0,0);CHKERRQ(ierr); 934c216dbf3SHong Zhang } 935c216dbf3SHong Zhang ierr = PetscLogEventBegin(MAT_TransposeMatMultNumeric,P,A,0,0);CHKERRQ(ierr); 936c216dbf3SHong Zhang ierr = (*(*C)->ops->mattransposemultnumeric)(P,A,*C);CHKERRQ(ierr); 937c216dbf3SHong Zhang ierr = PetscLogEventEnd(MAT_TransposeMatMultNumeric,P,A,0,0);CHKERRQ(ierr); 938ab1b013aSHong Zhang PetscFunctionReturn(0); 939ab1b013aSHong Zhang } 940ab1b013aSHong Zhang 941c216dbf3SHong Zhang /* This routine only works when scall=MAT_REUSE_MATRIX! */ 942c216dbf3SHong Zhang #undef __FUNCT__ 943c216dbf3SHong Zhang #define __FUNCT__ "MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult" 944c216dbf3SHong Zhang PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult(Mat P,Mat A,Mat C) 945c216dbf3SHong Zhang { 946c216dbf3SHong Zhang PetscErrorCode ierr; 9472bbb1c24SHong Zhang Mat_MPIAIJ *c=(Mat_MPIAIJ*)C->data; 9482bbb1c24SHong Zhang Mat_PtAPMPI *ptap= c->ptap; 9492bbb1c24SHong Zhang Mat Pt=ptap->Pt; 950c216dbf3SHong Zhang 951c216dbf3SHong Zhang PetscFunctionBegin; 952c216dbf3SHong Zhang ierr = MatTranspose(P,MAT_REUSE_MATRIX,&Pt);CHKERRQ(ierr); 953c216dbf3SHong Zhang ierr = MatMatMultNumeric(Pt,A,C);CHKERRQ(ierr); 954f96d37fcSHong Zhang PetscFunctionReturn(0); 955f96d37fcSHong Zhang } 956f96d37fcSHong Zhang 9576da69ca6SHong Zhang /* Non-scalable version, use dense axpy */ 958f96d37fcSHong Zhang #undef __FUNCT__ 9596da69ca6SHong Zhang #define __FUNCT__ "MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable" 9606da69ca6SHong Zhang PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable(Mat P,Mat A,Mat C) 961f96d37fcSHong Zhang { 962c5af039cSHong Zhang PetscErrorCode ierr; 963c5af039cSHong Zhang Mat_Merge_SeqsToMPI *merge; 964497f5370SHong Zhang Mat_MPIAIJ *p =(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data; 965c5af039cSHong Zhang Mat_SeqAIJ *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data; 966c5af039cSHong Zhang Mat_PtAPMPI *ptap; 967d6ab1dc0SHong Zhang PetscInt *adj,*aJ; 968497f5370SHong Zhang PetscInt i,j,k,anz,pnz,row,*cj; 969e745005bSHong Zhang MatScalar *ada,*aval,*ca,valtmp; 970c5af039cSHong Zhang PetscInt am =A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n; 971ce94432eSBarry Smith MPI_Comm comm; 972c5af039cSHong Zhang PetscMPIInt size,rank,taga,*len_s; 973c5af039cSHong Zhang PetscInt *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci; 974c5af039cSHong Zhang PetscInt **buf_ri,**buf_rj; 975c5af039cSHong Zhang PetscInt cnz=0,*bj_i,*bi,*bj,bnz,nextcj; /* bi,bj,ba: local array of C(mpi mat) */ 976c5af039cSHong Zhang MPI_Request *s_waits,*r_waits; 977c5af039cSHong Zhang MPI_Status *status; 978c5af039cSHong Zhang MatScalar **abuf_r,*ba_i,*pA,*coa,*ba; 979d6ab1dc0SHong Zhang PetscInt *ai,*aj,*coi,*coj; 980a2ea699eSBarry Smith PetscInt *poJ,*pdJ; 981c5af039cSHong Zhang Mat A_loc; 982c5af039cSHong Zhang Mat_SeqAIJ *a_loc; 983f96d37fcSHong Zhang 984f96d37fcSHong Zhang PetscFunctionBegin; 985ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr); 986c5af039cSHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 987c5af039cSHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 988c5af039cSHong Zhang 989c5af039cSHong Zhang ptap = c->ptap; 990c5af039cSHong Zhang merge = ptap->merge; 991c5af039cSHong Zhang 992c5af039cSHong Zhang /* 2) compute numeric C_seq = P_loc^T*A_loc*P - dominating part */ 993c5af039cSHong Zhang /*--------------------------------------------------------------*/ 994c5af039cSHong Zhang /* get data from symbolic products */ 995c5af039cSHong Zhang coi = merge->coi; coj = merge->coj; 996854ce69bSBarry Smith ierr = PetscCalloc1(coi[pon]+1,&coa);CHKERRQ(ierr); 997c5af039cSHong Zhang 998c5af039cSHong Zhang bi = merge->bi; bj = merge->bj; 999c5af039cSHong Zhang owners = merge->rowmap->range; 1000854ce69bSBarry Smith ierr = PetscCalloc1(bi[cm]+1,&ba);CHKERRQ(ierr); 1001c5af039cSHong Zhang 1002c5af039cSHong Zhang /* get A_loc by taking all local rows of A */ 1003c5af039cSHong Zhang A_loc = ptap->A_loc; 1004c5af039cSHong Zhang ierr = MatMPIAIJGetLocalMat(A,MAT_REUSE_MATRIX,&A_loc);CHKERRQ(ierr); 1005c5af039cSHong Zhang a_loc = (Mat_SeqAIJ*)(A_loc)->data; 1006d6ab1dc0SHong Zhang ai = a_loc->i; 1007d6ab1dc0SHong Zhang aj = a_loc->j; 1008c5af039cSHong Zhang 1009854ce69bSBarry Smith ierr = PetscCalloc1(A->cmap->N,&aval);CHKERRQ(ierr); /* non-scalable!!! */ 1010c5af039cSHong Zhang 1011c5af039cSHong Zhang for (i=0; i<am; i++) { 1012e745005bSHong Zhang /* 2-a) put A[i,:] to dense array aval */ 1013d6ab1dc0SHong Zhang anz = ai[i+1] - ai[i]; 1014d6ab1dc0SHong Zhang adj = aj + ai[i]; 1015d6ab1dc0SHong Zhang ada = a_loc->a + ai[i]; 1016c5af039cSHong Zhang for (j=0; j<anz; j++) { 1017e745005bSHong Zhang aval[adj[j]] = ada[j]; 1018c5af039cSHong Zhang } 1019c5af039cSHong Zhang 1020c5af039cSHong Zhang /* 2-b) Compute Cseq = P_loc[i,:]^T*A[i,:] using outer product */ 1021c5af039cSHong Zhang /*--------------------------------------------------------------*/ 1022c5af039cSHong Zhang /* put the value into Co=(p->B)^T*A (off-diagonal part, send to others) */ 1023c5af039cSHong Zhang pnz = po->i[i+1] - po->i[i]; 1024c5af039cSHong Zhang poJ = po->j + po->i[i]; 1025c5af039cSHong Zhang pA = po->a + po->i[i]; 1026c5af039cSHong Zhang for (j=0; j<pnz; j++) { 1027c5af039cSHong Zhang row = poJ[j]; 1028c5af039cSHong Zhang cnz = coi[row+1] - coi[row]; 1029c5af039cSHong Zhang cj = coj + coi[row]; 1030c5af039cSHong Zhang ca = coa + coi[row]; 1031c5af039cSHong Zhang /* perform dense axpy */ 1032c5af039cSHong Zhang valtmp = pA[j]; 1033c5af039cSHong Zhang for (k=0; k<cnz; k++) { 1034e745005bSHong Zhang ca[k] += valtmp*aval[cj[k]]; 1035c5af039cSHong Zhang } 1036c5af039cSHong Zhang ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr); 1037c5af039cSHong Zhang } 1038c5af039cSHong Zhang 1039c5af039cSHong Zhang /* put the value into Cd (diagonal part) */ 1040c5af039cSHong Zhang pnz = pd->i[i+1] - pd->i[i]; 1041c5af039cSHong Zhang pdJ = pd->j + pd->i[i]; 1042c5af039cSHong Zhang pA = pd->a + pd->i[i]; 1043c5af039cSHong Zhang for (j=0; j<pnz; j++) { 1044c5af039cSHong Zhang row = pdJ[j]; 1045c5af039cSHong Zhang cnz = bi[row+1] - bi[row]; 1046c5af039cSHong Zhang cj = bj + bi[row]; 1047c5af039cSHong Zhang ca = ba + bi[row]; 1048c5af039cSHong Zhang /* perform dense axpy */ 1049c5af039cSHong Zhang valtmp = pA[j]; 1050c5af039cSHong Zhang for (k=0; k<cnz; k++) { 1051e745005bSHong Zhang ca[k] += valtmp*aval[cj[k]]; 1052c5af039cSHong Zhang } 1053c5af039cSHong Zhang ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr); 1054c5af039cSHong Zhang } 1055c5af039cSHong Zhang 1056d6ab1dc0SHong Zhang /* zero the current row of Pt*A */ 1057d6ab1dc0SHong Zhang aJ = aj + ai[i]; 1058e745005bSHong Zhang for (k=0; k<anz; k++) aval[aJ[k]] = 0.0; 1059c5af039cSHong Zhang } 1060c5af039cSHong Zhang 1061c5af039cSHong Zhang /* 3) send and recv matrix values coa */ 1062c5af039cSHong Zhang /*------------------------------------*/ 1063c5af039cSHong Zhang buf_ri = merge->buf_ri; 1064c5af039cSHong Zhang buf_rj = merge->buf_rj; 1065c5af039cSHong Zhang len_s = merge->len_s; 1066c5af039cSHong Zhang ierr = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr); 1067c5af039cSHong Zhang ierr = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr); 1068c5af039cSHong Zhang 1069dcca6d9dSJed Brown ierr = PetscMalloc2(merge->nsend+1,&s_waits,size,&status);CHKERRQ(ierr); 1070c5af039cSHong Zhang for (proc=0,k=0; proc<size; proc++) { 1071c5af039cSHong Zhang if (!len_s[proc]) continue; 1072c5af039cSHong Zhang i = merge->owners_co[proc]; 1073c5af039cSHong Zhang ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr); 1074c5af039cSHong Zhang k++; 1075c5af039cSHong Zhang } 1076c5af039cSHong Zhang if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);} 1077c5af039cSHong Zhang if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);} 1078c5af039cSHong Zhang 1079c5af039cSHong Zhang ierr = PetscFree2(s_waits,status);CHKERRQ(ierr); 1080c5af039cSHong Zhang ierr = PetscFree(r_waits);CHKERRQ(ierr); 1081c5af039cSHong Zhang ierr = PetscFree(coa);CHKERRQ(ierr); 1082c5af039cSHong Zhang 1083c5af039cSHong Zhang /* 4) insert local Cseq and received values into Cmpi */ 1084c5af039cSHong Zhang /*----------------------------------------------------*/ 1085dcca6d9dSJed Brown ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr); 1086c5af039cSHong Zhang for (k=0; k<merge->nrecv; k++) { 1087c36aecf5SHong Zhang buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 1088c5af039cSHong Zhang nrows = *(buf_ri_k[k]); 1089c5af039cSHong Zhang nextrow[k] = buf_ri_k[k]+1; /* next row number of k-th recved i-structure */ 1090c5af039cSHong Zhang nextci[k] = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure */ 1091c5af039cSHong Zhang } 1092c5af039cSHong Zhang 1093c5af039cSHong Zhang for (i=0; i<cm; i++) { 1094c5af039cSHong Zhang row = owners[rank] + i; /* global row index of C_seq */ 1095c5af039cSHong Zhang bj_i = bj + bi[i]; /* col indices of the i-th row of C */ 1096c5af039cSHong Zhang ba_i = ba + bi[i]; 1097c5af039cSHong Zhang bnz = bi[i+1] - bi[i]; 1098c5af039cSHong Zhang /* add received vals into ba */ 1099c5af039cSHong Zhang for (k=0; k<merge->nrecv; k++) { /* k-th received message */ 1100c5af039cSHong Zhang /* i-th row */ 1101c5af039cSHong Zhang if (i == *nextrow[k]) { 1102c5af039cSHong Zhang cnz = *(nextci[k]+1) - *nextci[k]; 1103c5af039cSHong Zhang cj = buf_rj[k] + *(nextci[k]); 1104c5af039cSHong Zhang ca = abuf_r[k] + *(nextci[k]); 1105c5af039cSHong Zhang nextcj = 0; 1106c5af039cSHong Zhang for (j=0; nextcj<cnz; j++) { 1107c5af039cSHong Zhang if (bj_i[j] == cj[nextcj]) { /* bcol == ccol */ 1108c5af039cSHong Zhang ba_i[j] += ca[nextcj++]; 1109c5af039cSHong Zhang } 1110c5af039cSHong Zhang } 1111c5af039cSHong Zhang nextrow[k]++; nextci[k]++; 1112c5af039cSHong Zhang ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr); 1113c5af039cSHong Zhang } 1114c5af039cSHong Zhang } 1115c5af039cSHong Zhang ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr); 1116c5af039cSHong Zhang } 1117c5af039cSHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1118c5af039cSHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1119c5af039cSHong Zhang 1120c5af039cSHong Zhang ierr = PetscFree(ba);CHKERRQ(ierr); 1121c5af039cSHong Zhang ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr); 1122c5af039cSHong Zhang ierr = PetscFree(abuf_r);CHKERRQ(ierr); 1123c5af039cSHong Zhang ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr); 1124e745005bSHong Zhang ierr = PetscFree(aval);CHKERRQ(ierr); 1125f96d37fcSHong Zhang PetscFunctionReturn(0); 1126f96d37fcSHong Zhang } 1127f96d37fcSHong Zhang 1128c9ba551fSBarry Smith PetscErrorCode MatDuplicate_MPIAIJ_MatPtAP(Mat, MatDuplicateOption,Mat*); 1129c5af039cSHong Zhang /* This routine is modified from MatPtAPSymbolic_MPIAIJ_MPIAIJ() */ 1130f96d37fcSHong Zhang #undef __FUNCT__ 11312bbb1c24SHong Zhang #define __FUNCT__ "MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable" 11322bbb1c24SHong Zhang PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat P,Mat A,PetscReal fill,Mat *C) 1133f96d37fcSHong Zhang { 1134f96d37fcSHong Zhang PetscErrorCode ierr; 11354e938277SHong Zhang Mat Cmpi,A_loc,POt,PDt; 1136f96d37fcSHong Zhang Mat_PtAPMPI *ptap; 11370298fd71SBarry Smith PetscFreeSpaceList free_space=NULL,current_space=NULL; 1138f96d37fcSHong Zhang Mat_MPIAIJ *p =(Mat_MPIAIJ*)P->data,*c; 1139f96d37fcSHong Zhang PetscInt *pdti,*pdtj,*poti,*potj,*ptJ; 1140f96d37fcSHong Zhang PetscInt nnz; 11414e938277SHong Zhang PetscInt *lnk,*owners_co,*coi,*coj,i,k,pnz,row; 1142497f5370SHong Zhang PetscInt am=A->rmap->n,pn=P->cmap->n; 1143f96d37fcSHong Zhang PetscBT lnkbt; 1144ce94432eSBarry Smith MPI_Comm comm; 1145f96d37fcSHong Zhang PetscMPIInt size,rank,tagi,tagj,*len_si,*len_s,*len_ri; 1146f96d37fcSHong Zhang PetscInt **buf_rj,**buf_ri,**buf_ri_k; 1147f96d37fcSHong Zhang PetscInt len,proc,*dnz,*onz,*owners; 1148f96d37fcSHong Zhang PetscInt nzi,*bi,*bj; 1149f96d37fcSHong Zhang PetscInt nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci; 1150f96d37fcSHong Zhang MPI_Request *swaits,*rwaits; 1151f96d37fcSHong Zhang MPI_Status *sstatus,rstatus; 1152f96d37fcSHong Zhang Mat_Merge_SeqsToMPI *merge; 1153d6ab1dc0SHong Zhang PetscInt *ai,*aj,*Jptr,anz,*prmap=p->garray,pon,nspacedouble=0,j; 1154f96d37fcSHong Zhang PetscReal afill =1.0,afill_tmp; 1155c36aecf5SHong Zhang PetscInt rstart = P->cmap->rstart,rmax,aN=A->cmap->N,Crmax; 1156f96d37fcSHong Zhang PetscScalar *vals; 11574e938277SHong Zhang Mat_SeqAIJ *a_loc, *pdt,*pot; 1158f96d37fcSHong Zhang 1159f96d37fcSHong Zhang PetscFunctionBegin; 1160ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 1161c5af039cSHong Zhang /* check if matrix local sizes are compatible */ 1162c5af039cSHong Zhang if (A->rmap->rstart != P->rmap->rstart || A->rmap->rend != P->rmap->rend) { 1163c5af039cSHong 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); 1164c5af039cSHong Zhang } 1165c5af039cSHong Zhang 1166f96d37fcSHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1167f96d37fcSHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 1168f96d37fcSHong Zhang 1169f96d37fcSHong Zhang /* create struct Mat_PtAPMPI and attached it to C later */ 1170b00a9115SJed Brown ierr = PetscNew(&ptap);CHKERRQ(ierr); 1171f96d37fcSHong Zhang 1172f96d37fcSHong Zhang /* get A_loc by taking all local rows of A */ 1173f96d37fcSHong Zhang ierr = MatMPIAIJGetLocalMat(A,MAT_INITIAL_MATRIX,&A_loc);CHKERRQ(ierr); 11742205254eSKarl Rupp 1175c5af039cSHong Zhang ptap->A_loc = A_loc; 11762205254eSKarl Rupp 11771c7d5954SHong Zhang a_loc = (Mat_SeqAIJ*)(A_loc)->data; 1178d6ab1dc0SHong Zhang ai = a_loc->i; 1179d6ab1dc0SHong Zhang aj = a_loc->j; 1180f96d37fcSHong Zhang 1181f96d37fcSHong Zhang /* determine symbolic Co=(p->B)^T*A - send to others */ 1182f96d37fcSHong Zhang /*----------------------------------------------------*/ 11834e938277SHong Zhang ierr = MatTransposeSymbolic_SeqAIJ(p->A,&PDt);CHKERRQ(ierr); 11844e938277SHong Zhang pdt = (Mat_SeqAIJ*)PDt->data; 11854e938277SHong Zhang pdti = pdt->i; pdtj = pdt->j; 11864e938277SHong Zhang 11874e938277SHong Zhang ierr = MatTransposeSymbolic_SeqAIJ(p->B,&POt);CHKERRQ(ierr); 11884e938277SHong Zhang pot = (Mat_SeqAIJ*)POt->data; 11894e938277SHong Zhang poti = pot->i; potj = pot->j; 1190f96d37fcSHong Zhang 1191f96d37fcSHong Zhang /* then, compute symbolic Co = (p->B)^T*A */ 11922205254eSKarl Rupp pon = (p->B)->cmap->n; /* total num of rows to be sent to other processors >= (num of nonzero rows of C_seq) - pn */ 1193854ce69bSBarry Smith ierr = PetscMalloc1(pon+1,&coi);CHKERRQ(ierr); 1194f96d37fcSHong Zhang coi[0] = 0; 1195f96d37fcSHong Zhang 1196f96d37fcSHong Zhang /* set initial free space to be fill*(nnz(p->B) + nnz(A)) */ 1197d6ab1dc0SHong Zhang nnz = fill*(poti[pon] + ai[am]); 1198a2ea699eSBarry Smith ierr = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr); 1199f96d37fcSHong Zhang current_space = free_space; 120019f0e57aSHong Zhang 1201f96d37fcSHong Zhang /* create and initialize a linked list */ 12024e938277SHong Zhang i = PetscMax(pdt->rmax,pot->rmax); 1203c36aecf5SHong Zhang Crmax = i*a_loc->rmax*size; 12044e938277SHong Zhang if (!Crmax || Crmax > aN) Crmax = aN; 12054e938277SHong Zhang ierr = PetscLLCondensedCreate(Crmax,aN,&lnk,&lnkbt);CHKERRQ(ierr); 1206f96d37fcSHong Zhang 1207f96d37fcSHong Zhang for (i=0; i<pon; i++) { 1208f96d37fcSHong Zhang pnz = poti[i+1] - poti[i]; 1209f96d37fcSHong Zhang ptJ = potj + poti[i]; 1210f96d37fcSHong Zhang for (j=0; j<pnz; j++) { 1211f96d37fcSHong Zhang row = ptJ[j]; /* row of A_loc == col of Pot */ 1212d6ab1dc0SHong Zhang anz = ai[row+1] - ai[row]; 1213d6ab1dc0SHong Zhang Jptr = aj + ai[row]; 1214f96d37fcSHong Zhang /* add non-zero cols of AP into the sorted linked list lnk */ 1215d6ab1dc0SHong Zhang ierr = PetscLLCondensedAddSorted(anz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 1216f96d37fcSHong Zhang } 12174e938277SHong Zhang nnz = lnk[0]; 1218f96d37fcSHong Zhang 1219f96d37fcSHong Zhang /* If free space is not available, double the total space in the list */ 1220f96d37fcSHong Zhang if (current_space->local_remaining<nnz) { 1221f96d37fcSHong Zhang ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 1222f96d37fcSHong Zhang nspacedouble++; 1223f96d37fcSHong Zhang } 1224f96d37fcSHong Zhang 1225f96d37fcSHong Zhang /* Copy data into free space, and zero out denserows */ 12264e938277SHong Zhang ierr = PetscLLCondensedClean(aN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 12272205254eSKarl Rupp 1228f96d37fcSHong Zhang current_space->array += nnz; 1229f96d37fcSHong Zhang current_space->local_used += nnz; 1230f96d37fcSHong Zhang current_space->local_remaining -= nnz; 12312205254eSKarl Rupp 1232f96d37fcSHong Zhang coi[i+1] = coi[i] + nnz; 1233f96d37fcSHong Zhang } 1234f96d37fcSHong Zhang 1235854ce69bSBarry Smith ierr = PetscMalloc1(coi[pon]+1,&coj);CHKERRQ(ierr); 1236f96d37fcSHong Zhang ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr); 12372205254eSKarl Rupp 1238118727c9SMark F. Adams afill_tmp = (PetscReal)coi[pon]/(poti[pon] + ai[am]+1); 1239f96d37fcSHong Zhang if (afill_tmp > afill) afill = afill_tmp; 1240f96d37fcSHong Zhang 1241f96d37fcSHong Zhang /* send j-array (coj) of Co to other processors */ 1242f96d37fcSHong Zhang /*----------------------------------------------*/ 1243f96d37fcSHong Zhang /* determine row ownership */ 1244b00a9115SJed Brown ierr = PetscNew(&merge);CHKERRQ(ierr); 1245f96d37fcSHong Zhang ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr); 12462205254eSKarl Rupp 1247f96d37fcSHong Zhang merge->rowmap->n = pn; 1248f96d37fcSHong Zhang merge->rowmap->bs = 1; 12492205254eSKarl Rupp 1250f96d37fcSHong Zhang ierr = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr); 1251f96d37fcSHong Zhang owners = merge->rowmap->range; 1252f96d37fcSHong Zhang 1253f96d37fcSHong Zhang /* determine the number of messages to send, their lengths */ 12541795a4d1SJed Brown ierr = PetscCalloc1(size,&len_si);CHKERRQ(ierr); 1255785e854fSJed Brown ierr = PetscMalloc1(size,&merge->len_s);CHKERRQ(ierr); 12562205254eSKarl Rupp 1257f96d37fcSHong Zhang len_s = merge->len_s; 1258f96d37fcSHong Zhang merge->nsend = 0; 1259f96d37fcSHong Zhang 1260854ce69bSBarry Smith ierr = PetscMalloc1(size+2,&owners_co);CHKERRQ(ierr); 1261f96d37fcSHong Zhang ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 1262f96d37fcSHong Zhang 1263f96d37fcSHong Zhang proc = 0; 1264f96d37fcSHong Zhang for (i=0; i<pon; i++) { 1265f96d37fcSHong Zhang while (prmap[i] >= owners[proc+1]) proc++; 1266f96d37fcSHong Zhang len_si[proc]++; /* num of rows in Co to be sent to [proc] */ 1267f96d37fcSHong Zhang len_s[proc] += coi[i+1] - coi[i]; 1268f96d37fcSHong Zhang } 1269f96d37fcSHong Zhang 1270f96d37fcSHong Zhang len = 0; /* max length of buf_si[] */ 1271f96d37fcSHong Zhang owners_co[0] = 0; 1272f96d37fcSHong Zhang for (proc=0; proc<size; proc++) { 1273f96d37fcSHong Zhang owners_co[proc+1] = owners_co[proc] + len_si[proc]; 1274f96d37fcSHong Zhang if (len_si[proc]) { 1275f96d37fcSHong Zhang merge->nsend++; 1276f96d37fcSHong Zhang len_si[proc] = 2*(len_si[proc] + 1); 1277f96d37fcSHong Zhang len += len_si[proc]; 1278f96d37fcSHong Zhang } 1279f96d37fcSHong Zhang } 1280f96d37fcSHong Zhang 1281f96d37fcSHong Zhang /* determine the number and length of messages to receive for coi and coj */ 12820298fd71SBarry Smith ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&merge->nrecv);CHKERRQ(ierr); 1283f96d37fcSHong Zhang ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr); 1284f96d37fcSHong Zhang 1285f96d37fcSHong Zhang /* post the Irecv and Isend of coj */ 1286f96d37fcSHong Zhang ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr); 1287f96d37fcSHong Zhang ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr); 1288854ce69bSBarry Smith ierr = PetscMalloc1(merge->nsend+1,&swaits);CHKERRQ(ierr); 1289f96d37fcSHong Zhang for (proc=0, k=0; proc<size; proc++) { 1290f96d37fcSHong Zhang if (!len_s[proc]) continue; 1291f96d37fcSHong Zhang i = owners_co[proc]; 1292f96d37fcSHong Zhang ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr); 1293f96d37fcSHong Zhang k++; 1294f96d37fcSHong Zhang } 1295f96d37fcSHong Zhang 1296f96d37fcSHong Zhang /* receives and sends of coj are complete */ 1297785e854fSJed Brown ierr = PetscMalloc1(size,&sstatus);CHKERRQ(ierr); 1298f96d37fcSHong Zhang for (i=0; i<merge->nrecv; i++) { 1299c280ed6aSJed Brown PetscMPIInt icompleted; 1300c280ed6aSJed Brown ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 1301f96d37fcSHong Zhang } 1302f96d37fcSHong Zhang ierr = PetscFree(rwaits);CHKERRQ(ierr); 1303f96d37fcSHong Zhang if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);} 1304f96d37fcSHong Zhang 1305f96d37fcSHong Zhang /* send and recv coi */ 1306f96d37fcSHong Zhang /*-------------------*/ 1307f96d37fcSHong Zhang ierr = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr); 1308f96d37fcSHong Zhang ierr = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr); 1309854ce69bSBarry Smith ierr = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr); 1310f96d37fcSHong Zhang buf_si = buf_s; /* points to the beginning of k-th msg to be sent */ 1311f96d37fcSHong Zhang for (proc=0,k=0; proc<size; proc++) { 1312f96d37fcSHong Zhang if (!len_s[proc]) continue; 1313f96d37fcSHong Zhang /* form outgoing message for i-structure: 1314f96d37fcSHong Zhang buf_si[0]: nrows to be sent 1315f96d37fcSHong Zhang [1:nrows]: row index (global) 1316f96d37fcSHong Zhang [nrows+1:2*nrows+1]: i-structure index 1317f96d37fcSHong Zhang */ 1318f96d37fcSHong Zhang /*-------------------------------------------*/ 1319f96d37fcSHong Zhang nrows = len_si[proc]/2 - 1; 1320f96d37fcSHong Zhang buf_si_i = buf_si + nrows+1; 1321f96d37fcSHong Zhang buf_si[0] = nrows; 1322f96d37fcSHong Zhang buf_si_i[0] = 0; 1323f96d37fcSHong Zhang nrows = 0; 1324f96d37fcSHong Zhang for (i=owners_co[proc]; i<owners_co[proc+1]; i++) { 1325f96d37fcSHong Zhang nzi = coi[i+1] - coi[i]; 1326f96d37fcSHong Zhang buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */ 1327f96d37fcSHong Zhang buf_si[nrows+1] = prmap[i] -owners[proc]; /* local row index */ 1328f96d37fcSHong Zhang nrows++; 1329f96d37fcSHong Zhang } 1330f96d37fcSHong Zhang ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr); 1331f96d37fcSHong Zhang k++; 1332f96d37fcSHong Zhang buf_si += len_si[proc]; 1333f96d37fcSHong Zhang } 1334f96d37fcSHong Zhang i = merge->nrecv; 1335f96d37fcSHong Zhang while (i--) { 1336c280ed6aSJed Brown PetscMPIInt icompleted; 1337c280ed6aSJed Brown ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 1338f96d37fcSHong Zhang } 1339f96d37fcSHong Zhang ierr = PetscFree(rwaits);CHKERRQ(ierr); 1340f96d37fcSHong Zhang if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);} 1341f96d37fcSHong Zhang ierr = PetscFree(len_si);CHKERRQ(ierr); 1342f96d37fcSHong Zhang ierr = PetscFree(len_ri);CHKERRQ(ierr); 1343f96d37fcSHong Zhang ierr = PetscFree(swaits);CHKERRQ(ierr); 1344f96d37fcSHong Zhang ierr = PetscFree(sstatus);CHKERRQ(ierr); 1345f96d37fcSHong Zhang ierr = PetscFree(buf_s);CHKERRQ(ierr); 1346f96d37fcSHong Zhang 1347f96d37fcSHong Zhang /* compute the local portion of C (mpi mat) */ 1348f96d37fcSHong Zhang /*------------------------------------------*/ 1349f96d37fcSHong Zhang /* allocate bi array and free space for accumulating nonzero column info */ 1350854ce69bSBarry Smith ierr = PetscMalloc1(pn+1,&bi);CHKERRQ(ierr); 1351f96d37fcSHong Zhang bi[0] = 0; 1352f96d37fcSHong Zhang 1353c36aecf5SHong Zhang /* set initial free space to be fill*(nnz(P) + nnz(A)) */ 1354d6ab1dc0SHong Zhang nnz = fill*(pdti[pn] + poti[pon] + ai[am]); 1355a2ea699eSBarry Smith ierr = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr); 1356f96d37fcSHong Zhang current_space = free_space; 1357f96d37fcSHong Zhang 1358dcca6d9dSJed Brown ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr); 1359f96d37fcSHong Zhang for (k=0; k<merge->nrecv; k++) { 1360f96d37fcSHong Zhang buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 1361f96d37fcSHong Zhang nrows = *buf_ri_k[k]; 1362f96d37fcSHong Zhang nextrow[k] = buf_ri_k[k] + 1; /* next row number of k-th recved i-structure */ 1363f96d37fcSHong Zhang nextci[k] = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure */ 1364f96d37fcSHong Zhang } 1365f96d37fcSHong Zhang 13661c7d5954SHong Zhang ierr = MatPreallocateInitialize(comm,pn,A->cmap->n,dnz,onz);CHKERRQ(ierr); 1367f96d37fcSHong Zhang rmax = 0; 1368f96d37fcSHong Zhang for (i=0; i<pn; i++) { 1369f96d37fcSHong Zhang /* add pdt[i,:]*AP into lnk */ 1370f96d37fcSHong Zhang pnz = pdti[i+1] - pdti[i]; 1371f96d37fcSHong Zhang ptJ = pdtj + pdti[i]; 1372f96d37fcSHong Zhang for (j=0; j<pnz; j++) { 1373f96d37fcSHong Zhang row = ptJ[j]; /* row of AP == col of Pt */ 1374d6ab1dc0SHong Zhang anz = ai[row+1] - ai[row]; 1375d6ab1dc0SHong Zhang Jptr = aj + ai[row]; 1376f96d37fcSHong Zhang /* add non-zero cols of AP into the sorted linked list lnk */ 1377d6ab1dc0SHong Zhang ierr = PetscLLCondensedAddSorted(anz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 1378f96d37fcSHong Zhang } 13794e938277SHong Zhang 1380f96d37fcSHong Zhang /* add received col data into lnk */ 1381f96d37fcSHong Zhang for (k=0; k<merge->nrecv; k++) { /* k-th received message */ 1382f96d37fcSHong Zhang if (i == *nextrow[k]) { /* i-th row */ 1383f96d37fcSHong Zhang nzi = *(nextci[k]+1) - *nextci[k]; 1384f96d37fcSHong Zhang Jptr = buf_rj[k] + *nextci[k]; 13854e938277SHong Zhang ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr); 1386f96d37fcSHong Zhang nextrow[k]++; nextci[k]++; 1387f96d37fcSHong Zhang } 1388f96d37fcSHong Zhang } 13894e938277SHong Zhang nnz = lnk[0]; 1390f96d37fcSHong Zhang 1391f96d37fcSHong Zhang /* if free space is not available, make more free space */ 1392f96d37fcSHong Zhang if (current_space->local_remaining<nnz) { 1393f96d37fcSHong Zhang ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 1394f96d37fcSHong Zhang nspacedouble++; 1395f96d37fcSHong Zhang } 1396f96d37fcSHong Zhang /* copy data into free space, then initialize lnk */ 13974e938277SHong Zhang ierr = PetscLLCondensedClean(aN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 1398f96d37fcSHong Zhang ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr); 13992205254eSKarl Rupp 1400f96d37fcSHong Zhang current_space->array += nnz; 1401f96d37fcSHong Zhang current_space->local_used += nnz; 1402f96d37fcSHong Zhang current_space->local_remaining -= nnz; 14032205254eSKarl Rupp 1404f96d37fcSHong Zhang bi[i+1] = bi[i] + nnz; 1405f96d37fcSHong Zhang if (nnz > rmax) rmax = nnz; 1406f96d37fcSHong Zhang } 1407f96d37fcSHong Zhang ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr); 1408f96d37fcSHong Zhang 1409854ce69bSBarry Smith ierr = PetscMalloc1(bi[pn]+1,&bj);CHKERRQ(ierr); 1410f96d37fcSHong Zhang ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr); 14112205254eSKarl Rupp 1412118727c9SMark F. Adams afill_tmp = (PetscReal)bi[pn]/(pdti[pn] + poti[pon] + ai[am]+1); 1413f96d37fcSHong Zhang if (afill_tmp > afill) afill = afill_tmp; 1414d6ab1dc0SHong Zhang ierr = PetscLLCondensedDestroy(lnk,lnkbt);CHKERRQ(ierr); 14154e938277SHong Zhang ierr = MatDestroy(&POt);CHKERRQ(ierr); 14164e938277SHong Zhang ierr = MatDestroy(&PDt);CHKERRQ(ierr); 1417f96d37fcSHong Zhang 14181c7d5954SHong Zhang /* create symbolic parallel matrix Cmpi - why cannot be assembled in Numeric part */ 14191c7d5954SHong Zhang /*----------------------------------------------------------------------------------*/ 14201795a4d1SJed Brown ierr = PetscCalloc1(rmax+1,&vals);CHKERRQ(ierr); 1421f96d37fcSHong Zhang 1422f96d37fcSHong Zhang ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr); 1423f96d37fcSHong Zhang ierr = MatSetSizes(Cmpi,pn,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 142433d57670SJed Brown ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(A->cmap->bs));CHKERRQ(ierr); 1425f96d37fcSHong Zhang ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr); 1426f96d37fcSHong Zhang ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr); 1427f96d37fcSHong Zhang ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 1428f96d37fcSHong Zhang ierr = MatSetBlockSize(Cmpi,1);CHKERRQ(ierr); 1429f96d37fcSHong Zhang for (i=0; i<pn; i++) { 1430f96d37fcSHong Zhang row = i + rstart; 1431f96d37fcSHong Zhang nnz = bi[i+1] - bi[i]; 1432f96d37fcSHong Zhang Jptr = bj + bi[i]; 1433f96d37fcSHong Zhang ierr = MatSetValues(Cmpi,1,&row,nnz,Jptr,vals,INSERT_VALUES);CHKERRQ(ierr); 1434f96d37fcSHong Zhang } 1435f96d37fcSHong Zhang ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1436f96d37fcSHong Zhang ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1437f96d37fcSHong Zhang ierr = PetscFree(vals);CHKERRQ(ierr); 1438f96d37fcSHong Zhang 1439f96d37fcSHong Zhang merge->bi = bi; 1440f96d37fcSHong Zhang merge->bj = bj; 1441f96d37fcSHong Zhang merge->coi = coi; 1442f96d37fcSHong Zhang merge->coj = coj; 1443f96d37fcSHong Zhang merge->buf_ri = buf_ri; 1444f96d37fcSHong Zhang merge->buf_rj = buf_rj; 1445f96d37fcSHong Zhang merge->owners_co = owners_co; 1446f96d37fcSHong Zhang merge->destroy = Cmpi->ops->destroy; 1447f96d37fcSHong Zhang merge->duplicate = Cmpi->ops->duplicate; 1448f96d37fcSHong Zhang 14496da69ca6SHong Zhang Cmpi->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable; 1450c5af039cSHong Zhang Cmpi->ops->destroy = MatDestroy_MPIAIJ_PtAP; 1451c9ba551fSBarry Smith Cmpi->ops->duplicate = MatDuplicate_MPIAIJ_MatPtAP; 1452f96d37fcSHong Zhang 1453f96d37fcSHong Zhang /* attach the supporting struct to Cmpi for reuse */ 1454f96d37fcSHong Zhang c = (Mat_MPIAIJ*)Cmpi->data; 1455f96d37fcSHong Zhang c->ptap = ptap; 14560298fd71SBarry Smith ptap->api = NULL; 14570298fd71SBarry Smith ptap->apj = NULL; 1458f96d37fcSHong Zhang ptap->merge = merge; 1459d6ab1dc0SHong Zhang ptap->rmax = rmax; 1460d6ab1dc0SHong Zhang 1461d6ab1dc0SHong Zhang *C = Cmpi; 1462d6ab1dc0SHong Zhang #if defined(PETSC_USE_INFO) 1463d6ab1dc0SHong Zhang if (bi[pn] != 0) { 146457622a8eSBarry Smith ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr); 146557622a8eSBarry Smith ierr = PetscInfo1(Cmpi,"Use MatTransposeMatMult(A,B,MatReuse,%g,&C) for best performance.\n",(double)afill);CHKERRQ(ierr); 1466d6ab1dc0SHong Zhang } else { 1467d6ab1dc0SHong Zhang ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr); 1468d6ab1dc0SHong Zhang } 1469d6ab1dc0SHong Zhang #endif 1470d6ab1dc0SHong Zhang PetscFunctionReturn(0); 1471d6ab1dc0SHong Zhang } 1472d6ab1dc0SHong Zhang 1473d6ab1dc0SHong Zhang #undef __FUNCT__ 14746da69ca6SHong Zhang #define __FUNCT__ "MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ" 14756da69ca6SHong Zhang PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ(Mat P,Mat A,Mat C) 1476d6ab1dc0SHong Zhang { 1477d6ab1dc0SHong Zhang PetscErrorCode ierr; 1478d6ab1dc0SHong Zhang Mat_Merge_SeqsToMPI *merge; 1479d6ab1dc0SHong Zhang Mat_MPIAIJ *p =(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data; 1480d6ab1dc0SHong Zhang Mat_SeqAIJ *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data; 1481d6ab1dc0SHong Zhang Mat_PtAPMPI *ptap; 1482e745005bSHong Zhang PetscInt *adj; 1483e745005bSHong Zhang PetscInt i,j,k,anz,pnz,row,*cj,nexta; 1484e745005bSHong Zhang MatScalar *ada,*ca,valtmp; 1485d6ab1dc0SHong Zhang PetscInt am =A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n; 1486ce94432eSBarry Smith MPI_Comm comm; 1487d6ab1dc0SHong Zhang PetscMPIInt size,rank,taga,*len_s; 1488d6ab1dc0SHong Zhang PetscInt *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci; 1489d6ab1dc0SHong Zhang PetscInt **buf_ri,**buf_rj; 1490d6ab1dc0SHong Zhang PetscInt cnz=0,*bj_i,*bi,*bj,bnz,nextcj; /* bi,bj,ba: local array of C(mpi mat) */ 1491d6ab1dc0SHong Zhang MPI_Request *s_waits,*r_waits; 1492d6ab1dc0SHong Zhang MPI_Status *status; 1493d6ab1dc0SHong Zhang MatScalar **abuf_r,*ba_i,*pA,*coa,*ba; 1494d6ab1dc0SHong Zhang PetscInt *ai,*aj,*coi,*coj; 1495a2ea699eSBarry Smith PetscInt *poJ,*pdJ; 1496d6ab1dc0SHong Zhang Mat A_loc; 1497d6ab1dc0SHong Zhang Mat_SeqAIJ *a_loc; 1498d6ab1dc0SHong Zhang 1499d6ab1dc0SHong Zhang PetscFunctionBegin; 1500ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr); 1501d6ab1dc0SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1502d6ab1dc0SHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 1503d6ab1dc0SHong Zhang 1504d6ab1dc0SHong Zhang ptap = c->ptap; 1505d6ab1dc0SHong Zhang merge = ptap->merge; 1506d6ab1dc0SHong Zhang 1507e745005bSHong Zhang /* 2) compute numeric C_seq = P_loc^T*A_loc */ 1508e745005bSHong Zhang /*------------------------------------------*/ 1509d6ab1dc0SHong Zhang /* get data from symbolic products */ 1510d6ab1dc0SHong Zhang coi = merge->coi; coj = merge->coj; 1511854ce69bSBarry Smith ierr = PetscCalloc1(coi[pon]+1,&coa);CHKERRQ(ierr); 1512d6ab1dc0SHong Zhang bi = merge->bi; bj = merge->bj; 1513d6ab1dc0SHong Zhang owners = merge->rowmap->range; 15141795a4d1SJed Brown ierr = PetscCalloc1(bi[cm]+1,&ba);CHKERRQ(ierr); 1515d6ab1dc0SHong Zhang 1516d6ab1dc0SHong Zhang /* get A_loc by taking all local rows of A */ 1517d6ab1dc0SHong Zhang A_loc = ptap->A_loc; 1518d6ab1dc0SHong Zhang ierr = MatMPIAIJGetLocalMat(A,MAT_REUSE_MATRIX,&A_loc);CHKERRQ(ierr); 1519d6ab1dc0SHong Zhang a_loc = (Mat_SeqAIJ*)(A_loc)->data; 1520d6ab1dc0SHong Zhang ai = a_loc->i; 1521d6ab1dc0SHong Zhang aj = a_loc->j; 1522d6ab1dc0SHong Zhang 1523d6ab1dc0SHong Zhang for (i=0; i<am; i++) { 1524d6ab1dc0SHong Zhang anz = ai[i+1] - ai[i]; 1525d6ab1dc0SHong Zhang adj = aj + ai[i]; 1526d6ab1dc0SHong Zhang ada = a_loc->a + ai[i]; 1527d6ab1dc0SHong Zhang 1528d6ab1dc0SHong Zhang /* 2-b) Compute Cseq = P_loc[i,:]^T*A[i,:] using outer product */ 1529e745005bSHong Zhang /*-------------------------------------------------------------*/ 1530d6ab1dc0SHong Zhang /* put the value into Co=(p->B)^T*A (off-diagonal part, send to others) */ 1531d6ab1dc0SHong Zhang pnz = po->i[i+1] - po->i[i]; 1532d6ab1dc0SHong Zhang poJ = po->j + po->i[i]; 1533d6ab1dc0SHong Zhang pA = po->a + po->i[i]; 1534d6ab1dc0SHong Zhang for (j=0; j<pnz; j++) { 1535d6ab1dc0SHong Zhang row = poJ[j]; 1536d6ab1dc0SHong Zhang cj = coj + coi[row]; 1537d6ab1dc0SHong Zhang ca = coa + coi[row]; 1538e745005bSHong Zhang /* perform sparse axpy */ 1539e745005bSHong Zhang nexta = 0; 1540d6ab1dc0SHong Zhang valtmp = pA[j]; 1541e745005bSHong Zhang for (k=0; nexta<anz; k++) { 1542e745005bSHong Zhang if (cj[k] == adj[nexta]) { 1543e745005bSHong Zhang ca[k] += valtmp*ada[nexta]; 1544e745005bSHong Zhang nexta++; 1545d6ab1dc0SHong Zhang } 1546e745005bSHong Zhang } 1547e745005bSHong Zhang ierr = PetscLogFlops(2.0*anz);CHKERRQ(ierr); 1548d6ab1dc0SHong Zhang } 1549d6ab1dc0SHong Zhang 1550d6ab1dc0SHong Zhang /* put the value into Cd (diagonal part) */ 1551d6ab1dc0SHong Zhang pnz = pd->i[i+1] - pd->i[i]; 1552d6ab1dc0SHong Zhang pdJ = pd->j + pd->i[i]; 1553d6ab1dc0SHong Zhang pA = pd->a + pd->i[i]; 1554d6ab1dc0SHong Zhang for (j=0; j<pnz; j++) { 1555d6ab1dc0SHong Zhang row = pdJ[j]; 1556d6ab1dc0SHong Zhang cj = bj + bi[row]; 1557d6ab1dc0SHong Zhang ca = ba + bi[row]; 1558e745005bSHong Zhang /* perform sparse axpy */ 1559e745005bSHong Zhang nexta = 0; 1560d6ab1dc0SHong Zhang valtmp = pA[j]; 1561e745005bSHong Zhang for (k=0; nexta<anz; k++) { 1562e745005bSHong Zhang if (cj[k] == adj[nexta]) { 1563e745005bSHong Zhang ca[k] += valtmp*ada[nexta]; 1564e745005bSHong Zhang nexta++; 1565d6ab1dc0SHong Zhang } 1566e745005bSHong Zhang } 1567e745005bSHong Zhang ierr = PetscLogFlops(2.0*anz);CHKERRQ(ierr); 1568d6ab1dc0SHong Zhang } 1569d6ab1dc0SHong Zhang } 1570d6ab1dc0SHong Zhang 1571d6ab1dc0SHong Zhang /* 3) send and recv matrix values coa */ 1572d6ab1dc0SHong Zhang /*------------------------------------*/ 1573d6ab1dc0SHong Zhang buf_ri = merge->buf_ri; 1574d6ab1dc0SHong Zhang buf_rj = merge->buf_rj; 1575d6ab1dc0SHong Zhang len_s = merge->len_s; 1576d6ab1dc0SHong Zhang ierr = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr); 1577d6ab1dc0SHong Zhang ierr = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr); 1578d6ab1dc0SHong Zhang 1579dcca6d9dSJed Brown ierr = PetscMalloc2(merge->nsend+1,&s_waits,size,&status);CHKERRQ(ierr); 1580d6ab1dc0SHong Zhang for (proc=0,k=0; proc<size; proc++) { 1581d6ab1dc0SHong Zhang if (!len_s[proc]) continue; 1582d6ab1dc0SHong Zhang i = merge->owners_co[proc]; 1583d6ab1dc0SHong Zhang ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr); 1584d6ab1dc0SHong Zhang k++; 1585d6ab1dc0SHong Zhang } 1586d6ab1dc0SHong Zhang if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);} 1587d6ab1dc0SHong Zhang if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);} 1588d6ab1dc0SHong Zhang 1589d6ab1dc0SHong Zhang ierr = PetscFree2(s_waits,status);CHKERRQ(ierr); 1590d6ab1dc0SHong Zhang ierr = PetscFree(r_waits);CHKERRQ(ierr); 1591d6ab1dc0SHong Zhang ierr = PetscFree(coa);CHKERRQ(ierr); 1592d6ab1dc0SHong Zhang 1593d6ab1dc0SHong Zhang /* 4) insert local Cseq and received values into Cmpi */ 1594d6ab1dc0SHong Zhang /*----------------------------------------------------*/ 1595dcca6d9dSJed Brown ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr); 1596d6ab1dc0SHong Zhang for (k=0; k<merge->nrecv; k++) { 1597e745005bSHong Zhang buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 1598d6ab1dc0SHong Zhang nrows = *(buf_ri_k[k]); 1599d6ab1dc0SHong Zhang nextrow[k] = buf_ri_k[k]+1; /* next row number of k-th recved i-structure */ 1600d6ab1dc0SHong Zhang nextci[k] = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure */ 1601d6ab1dc0SHong Zhang } 1602d6ab1dc0SHong Zhang 1603d6ab1dc0SHong Zhang for (i=0; i<cm; i++) { 1604d6ab1dc0SHong Zhang row = owners[rank] + i; /* global row index of C_seq */ 1605d6ab1dc0SHong Zhang bj_i = bj + bi[i]; /* col indices of the i-th row of C */ 1606d6ab1dc0SHong Zhang ba_i = ba + bi[i]; 1607d6ab1dc0SHong Zhang bnz = bi[i+1] - bi[i]; 1608d6ab1dc0SHong Zhang /* add received vals into ba */ 1609d6ab1dc0SHong Zhang for (k=0; k<merge->nrecv; k++) { /* k-th received message */ 1610d6ab1dc0SHong Zhang /* i-th row */ 1611d6ab1dc0SHong Zhang if (i == *nextrow[k]) { 1612d6ab1dc0SHong Zhang cnz = *(nextci[k]+1) - *nextci[k]; 1613d6ab1dc0SHong Zhang cj = buf_rj[k] + *(nextci[k]); 1614d6ab1dc0SHong Zhang ca = abuf_r[k] + *(nextci[k]); 1615d6ab1dc0SHong Zhang nextcj = 0; 1616d6ab1dc0SHong Zhang for (j=0; nextcj<cnz; j++) { 1617d6ab1dc0SHong Zhang if (bj_i[j] == cj[nextcj]) { /* bcol == ccol */ 1618d6ab1dc0SHong Zhang ba_i[j] += ca[nextcj++]; 1619d6ab1dc0SHong Zhang } 1620d6ab1dc0SHong Zhang } 1621d6ab1dc0SHong Zhang nextrow[k]++; nextci[k]++; 1622d6ab1dc0SHong Zhang ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr); 1623d6ab1dc0SHong Zhang } 1624d6ab1dc0SHong Zhang } 1625d6ab1dc0SHong Zhang ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr); 1626d6ab1dc0SHong Zhang } 1627d6ab1dc0SHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1628d6ab1dc0SHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1629d6ab1dc0SHong Zhang 1630d6ab1dc0SHong Zhang ierr = PetscFree(ba);CHKERRQ(ierr); 1631d6ab1dc0SHong Zhang ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr); 1632d6ab1dc0SHong Zhang ierr = PetscFree(abuf_r);CHKERRQ(ierr); 1633d6ab1dc0SHong Zhang ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr); 1634d6ab1dc0SHong Zhang PetscFunctionReturn(0); 1635d6ab1dc0SHong Zhang } 1636d6ab1dc0SHong Zhang 1637c9ba551fSBarry Smith PetscErrorCode MatDuplicate_MPIAIJ_MatPtAP(Mat, MatDuplicateOption,Mat*); 16382bbb1c24SHong Zhang /* This routine is modified from MatPtAPSymbolic_MPIAIJ_MPIAIJ(); 16392bbb1c24SHong Zhang differ from MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable in using LLCondensedCreate_Scalable() */ 1640d6ab1dc0SHong Zhang #undef __FUNCT__ 16416da69ca6SHong Zhang #define __FUNCT__ "MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ" 16426da69ca6SHong Zhang PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(Mat P,Mat A,PetscReal fill,Mat *C) 1643d6ab1dc0SHong Zhang { 1644d6ab1dc0SHong Zhang PetscErrorCode ierr; 1645d6ab1dc0SHong Zhang Mat Cmpi,A_loc,POt,PDt; 1646d6ab1dc0SHong Zhang Mat_PtAPMPI *ptap; 16470298fd71SBarry Smith PetscFreeSpaceList free_space=NULL,current_space=NULL; 1648d6ab1dc0SHong Zhang Mat_MPIAIJ *p =(Mat_MPIAIJ*)P->data,*c; 1649d6ab1dc0SHong Zhang PetscInt *pdti,*pdtj,*poti,*potj,*ptJ; 1650d6ab1dc0SHong Zhang PetscInt nnz; 1651d6ab1dc0SHong Zhang PetscInt *lnk,*owners_co,*coi,*coj,i,k,pnz,row; 1652d6ab1dc0SHong Zhang PetscInt am =A->rmap->n,pn=P->cmap->n; 1653ce94432eSBarry Smith MPI_Comm comm; 1654d6ab1dc0SHong Zhang PetscMPIInt size,rank,tagi,tagj,*len_si,*len_s,*len_ri; 1655d6ab1dc0SHong Zhang PetscInt **buf_rj,**buf_ri,**buf_ri_k; 1656d6ab1dc0SHong Zhang PetscInt len,proc,*dnz,*onz,*owners; 1657d6ab1dc0SHong Zhang PetscInt nzi,*bi,*bj; 1658d6ab1dc0SHong Zhang PetscInt nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci; 1659d6ab1dc0SHong Zhang MPI_Request *swaits,*rwaits; 1660d6ab1dc0SHong Zhang MPI_Status *sstatus,rstatus; 1661d6ab1dc0SHong Zhang Mat_Merge_SeqsToMPI *merge; 1662d6ab1dc0SHong Zhang PetscInt *ai,*aj,*Jptr,anz,*prmap=p->garray,pon,nspacedouble=0,j; 1663d6ab1dc0SHong Zhang PetscReal afill =1.0,afill_tmp; 1664c36aecf5SHong Zhang PetscInt rstart = P->cmap->rstart,rmax,aN=A->cmap->N,Crmax; 1665d6ab1dc0SHong Zhang PetscScalar *vals; 1666d6ab1dc0SHong Zhang Mat_SeqAIJ *a_loc, *pdt,*pot; 1667d6ab1dc0SHong Zhang 1668d6ab1dc0SHong Zhang PetscFunctionBegin; 1669ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 1670d6ab1dc0SHong Zhang /* check if matrix local sizes are compatible */ 1671d6ab1dc0SHong Zhang if (A->rmap->rstart != P->rmap->rstart || A->rmap->rend != P->rmap->rend) { 1672d6ab1dc0SHong 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); 1673d6ab1dc0SHong Zhang } 1674d6ab1dc0SHong Zhang 1675d6ab1dc0SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1676d6ab1dc0SHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 1677d6ab1dc0SHong Zhang 1678d6ab1dc0SHong Zhang /* create struct Mat_PtAPMPI and attached it to C later */ 1679b00a9115SJed Brown ierr = PetscNew(&ptap);CHKERRQ(ierr); 1680d6ab1dc0SHong Zhang 1681d6ab1dc0SHong Zhang /* get A_loc by taking all local rows of A */ 1682d6ab1dc0SHong Zhang ierr = MatMPIAIJGetLocalMat(A,MAT_INITIAL_MATRIX,&A_loc);CHKERRQ(ierr); 16832205254eSKarl Rupp 1684d6ab1dc0SHong Zhang ptap->A_loc = A_loc; 1685d6ab1dc0SHong Zhang a_loc = (Mat_SeqAIJ*)(A_loc)->data; 1686d6ab1dc0SHong Zhang ai = a_loc->i; 1687d6ab1dc0SHong Zhang aj = a_loc->j; 1688d6ab1dc0SHong Zhang 1689d6ab1dc0SHong Zhang /* determine symbolic Co=(p->B)^T*A - send to others */ 1690d6ab1dc0SHong Zhang /*----------------------------------------------------*/ 1691d6ab1dc0SHong Zhang ierr = MatTransposeSymbolic_SeqAIJ(p->A,&PDt);CHKERRQ(ierr); 1692d6ab1dc0SHong Zhang pdt = (Mat_SeqAIJ*)PDt->data; 1693d6ab1dc0SHong Zhang pdti = pdt->i; pdtj = pdt->j; 1694d6ab1dc0SHong Zhang 1695d6ab1dc0SHong Zhang ierr = MatTransposeSymbolic_SeqAIJ(p->B,&POt);CHKERRQ(ierr); 1696d6ab1dc0SHong Zhang pot = (Mat_SeqAIJ*)POt->data; 1697d6ab1dc0SHong Zhang poti = pot->i; potj = pot->j; 1698d6ab1dc0SHong Zhang 1699d6ab1dc0SHong Zhang /* then, compute symbolic Co = (p->B)^T*A */ 1700d6ab1dc0SHong Zhang pon = (p->B)->cmap->n; /* total num of rows to be sent to other processors 1701d6ab1dc0SHong Zhang >= (num of nonzero rows of C_seq) - pn */ 1702854ce69bSBarry Smith ierr = PetscMalloc1(pon+1,&coi);CHKERRQ(ierr); 1703d6ab1dc0SHong Zhang coi[0] = 0; 1704d6ab1dc0SHong Zhang 1705d6ab1dc0SHong Zhang /* set initial free space to be fill*(nnz(p->B) + nnz(A)) */ 1706d6ab1dc0SHong Zhang nnz = fill*(poti[pon] + ai[am]); 1707a2ea699eSBarry Smith ierr = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr); 1708d6ab1dc0SHong Zhang current_space = free_space; 170919f0e57aSHong Zhang 1710d6ab1dc0SHong Zhang /* create and initialize a linked list */ 1711d6ab1dc0SHong Zhang i = PetscMax(pdt->rmax,pot->rmax); 1712c36aecf5SHong Zhang Crmax = i*a_loc->rmax*size; /* non-scalable! */ 1713d6ab1dc0SHong Zhang if (!Crmax || Crmax > aN) Crmax = aN; 1714d6ab1dc0SHong Zhang ierr = PetscLLCondensedCreate_Scalable(Crmax,&lnk);CHKERRQ(ierr); 1715d6ab1dc0SHong Zhang 1716d6ab1dc0SHong Zhang for (i=0; i<pon; i++) { 1717d6ab1dc0SHong Zhang pnz = poti[i+1] - poti[i]; 1718d6ab1dc0SHong Zhang ptJ = potj + poti[i]; 1719d6ab1dc0SHong Zhang for (j=0; j<pnz; j++) { 1720d6ab1dc0SHong Zhang row = ptJ[j]; /* row of A_loc == col of Pot */ 1721d6ab1dc0SHong Zhang anz = ai[row+1] - ai[row]; 1722d6ab1dc0SHong Zhang Jptr = aj + ai[row]; 1723d6ab1dc0SHong Zhang /* add non-zero cols of AP into the sorted linked list lnk */ 1724d6ab1dc0SHong Zhang ierr = PetscLLCondensedAddSorted_Scalable(anz,Jptr,lnk);CHKERRQ(ierr); 1725d6ab1dc0SHong Zhang } 1726d6ab1dc0SHong Zhang nnz = lnk[0]; 1727d6ab1dc0SHong Zhang 1728d6ab1dc0SHong Zhang /* If free space is not available, double the total space in the list */ 1729d6ab1dc0SHong Zhang if (current_space->local_remaining<nnz) { 1730d6ab1dc0SHong Zhang ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 1731d6ab1dc0SHong Zhang nspacedouble++; 1732d6ab1dc0SHong Zhang } 1733d6ab1dc0SHong Zhang 1734d6ab1dc0SHong Zhang /* Copy data into free space, and zero out denserows */ 1735d6ab1dc0SHong Zhang ierr = PetscLLCondensedClean_Scalable(nnz,current_space->array,lnk);CHKERRQ(ierr); 17362205254eSKarl Rupp 1737d6ab1dc0SHong Zhang current_space->array += nnz; 1738d6ab1dc0SHong Zhang current_space->local_used += nnz; 1739d6ab1dc0SHong Zhang current_space->local_remaining -= nnz; 17402205254eSKarl Rupp 1741d6ab1dc0SHong Zhang coi[i+1] = coi[i] + nnz; 1742d6ab1dc0SHong Zhang } 1743d6ab1dc0SHong Zhang 1744854ce69bSBarry Smith ierr = PetscMalloc1(coi[pon]+1,&coj);CHKERRQ(ierr); 1745d6ab1dc0SHong Zhang ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr); 17462205254eSKarl Rupp 1747118727c9SMark F. Adams afill_tmp = (PetscReal)coi[pon]/(poti[pon] + ai[am]+1); 1748d6ab1dc0SHong Zhang if (afill_tmp > afill) afill = afill_tmp; 1749d6ab1dc0SHong Zhang 1750d6ab1dc0SHong Zhang /* send j-array (coj) of Co to other processors */ 1751d6ab1dc0SHong Zhang /*----------------------------------------------*/ 1752d6ab1dc0SHong Zhang /* determine row ownership */ 1753b00a9115SJed Brown ierr = PetscNew(&merge);CHKERRQ(ierr); 1754d6ab1dc0SHong Zhang ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr); 17552205254eSKarl Rupp 1756d6ab1dc0SHong Zhang merge->rowmap->n = pn; 1757d6ab1dc0SHong Zhang merge->rowmap->bs = 1; 17582205254eSKarl Rupp 1759d6ab1dc0SHong Zhang ierr = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr); 1760d6ab1dc0SHong Zhang owners = merge->rowmap->range; 1761d6ab1dc0SHong Zhang 1762d6ab1dc0SHong Zhang /* determine the number of messages to send, their lengths */ 17631795a4d1SJed Brown ierr = PetscCalloc1(size,&len_si);CHKERRQ(ierr); 1764785e854fSJed Brown ierr = PetscMalloc1(size,&merge->len_s);CHKERRQ(ierr); 17652205254eSKarl Rupp 1766d6ab1dc0SHong Zhang len_s = merge->len_s; 1767d6ab1dc0SHong Zhang merge->nsend = 0; 1768d6ab1dc0SHong Zhang 1769854ce69bSBarry Smith ierr = PetscMalloc1(size+2,&owners_co);CHKERRQ(ierr); 1770d6ab1dc0SHong Zhang ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 1771d6ab1dc0SHong Zhang 1772d6ab1dc0SHong Zhang proc = 0; 1773d6ab1dc0SHong Zhang for (i=0; i<pon; i++) { 1774d6ab1dc0SHong Zhang while (prmap[i] >= owners[proc+1]) proc++; 1775d6ab1dc0SHong Zhang len_si[proc]++; /* num of rows in Co to be sent to [proc] */ 1776d6ab1dc0SHong Zhang len_s[proc] += coi[i+1] - coi[i]; 1777d6ab1dc0SHong Zhang } 1778d6ab1dc0SHong Zhang 1779d6ab1dc0SHong Zhang len = 0; /* max length of buf_si[] */ 1780d6ab1dc0SHong Zhang owners_co[0] = 0; 1781d6ab1dc0SHong Zhang for (proc=0; proc<size; proc++) { 1782d6ab1dc0SHong Zhang owners_co[proc+1] = owners_co[proc] + len_si[proc]; 1783d6ab1dc0SHong Zhang if (len_si[proc]) { 1784d6ab1dc0SHong Zhang merge->nsend++; 1785d6ab1dc0SHong Zhang len_si[proc] = 2*(len_si[proc] + 1); 1786d6ab1dc0SHong Zhang len += len_si[proc]; 1787d6ab1dc0SHong Zhang } 1788d6ab1dc0SHong Zhang } 1789d6ab1dc0SHong Zhang 1790d6ab1dc0SHong Zhang /* determine the number and length of messages to receive for coi and coj */ 17910298fd71SBarry Smith ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&merge->nrecv);CHKERRQ(ierr); 1792d6ab1dc0SHong Zhang ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr); 1793d6ab1dc0SHong Zhang 1794d6ab1dc0SHong Zhang /* post the Irecv and Isend of coj */ 1795d6ab1dc0SHong Zhang ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr); 1796d6ab1dc0SHong Zhang ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr); 1797854ce69bSBarry Smith ierr = PetscMalloc1(merge->nsend+1,&swaits);CHKERRQ(ierr); 1798d6ab1dc0SHong Zhang for (proc=0, k=0; proc<size; proc++) { 1799d6ab1dc0SHong Zhang if (!len_s[proc]) continue; 1800d6ab1dc0SHong Zhang i = owners_co[proc]; 1801d6ab1dc0SHong Zhang ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr); 1802d6ab1dc0SHong Zhang k++; 1803d6ab1dc0SHong Zhang } 1804d6ab1dc0SHong Zhang 1805d6ab1dc0SHong Zhang /* receives and sends of coj are complete */ 1806785e854fSJed Brown ierr = PetscMalloc1(size,&sstatus);CHKERRQ(ierr); 1807d6ab1dc0SHong Zhang for (i=0; i<merge->nrecv; i++) { 1808c280ed6aSJed Brown PetscMPIInt icompleted; 1809c280ed6aSJed Brown ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 1810d6ab1dc0SHong Zhang } 1811d6ab1dc0SHong Zhang ierr = PetscFree(rwaits);CHKERRQ(ierr); 1812d6ab1dc0SHong Zhang if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);} 1813d6ab1dc0SHong Zhang 1814d6ab1dc0SHong Zhang /* send and recv coi */ 1815d6ab1dc0SHong Zhang /*-------------------*/ 1816d6ab1dc0SHong Zhang ierr = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr); 1817d6ab1dc0SHong Zhang ierr = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr); 1818854ce69bSBarry Smith ierr = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr); 1819d6ab1dc0SHong Zhang buf_si = buf_s; /* points to the beginning of k-th msg to be sent */ 1820d6ab1dc0SHong Zhang for (proc=0,k=0; proc<size; proc++) { 1821d6ab1dc0SHong Zhang if (!len_s[proc]) continue; 1822d6ab1dc0SHong Zhang /* form outgoing message for i-structure: 1823d6ab1dc0SHong Zhang buf_si[0]: nrows to be sent 1824d6ab1dc0SHong Zhang [1:nrows]: row index (global) 1825d6ab1dc0SHong Zhang [nrows+1:2*nrows+1]: i-structure index 1826d6ab1dc0SHong Zhang */ 1827d6ab1dc0SHong Zhang /*-------------------------------------------*/ 1828d6ab1dc0SHong Zhang nrows = len_si[proc]/2 - 1; 1829d6ab1dc0SHong Zhang buf_si_i = buf_si + nrows+1; 1830d6ab1dc0SHong Zhang buf_si[0] = nrows; 1831d6ab1dc0SHong Zhang buf_si_i[0] = 0; 1832d6ab1dc0SHong Zhang nrows = 0; 1833d6ab1dc0SHong Zhang for (i=owners_co[proc]; i<owners_co[proc+1]; i++) { 1834d6ab1dc0SHong Zhang nzi = coi[i+1] - coi[i]; 1835d6ab1dc0SHong Zhang buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */ 1836d6ab1dc0SHong Zhang buf_si[nrows+1] = prmap[i] -owners[proc]; /* local row index */ 1837d6ab1dc0SHong Zhang nrows++; 1838d6ab1dc0SHong Zhang } 1839d6ab1dc0SHong Zhang ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr); 1840d6ab1dc0SHong Zhang k++; 1841d6ab1dc0SHong Zhang buf_si += len_si[proc]; 1842d6ab1dc0SHong Zhang } 1843d6ab1dc0SHong Zhang i = merge->nrecv; 1844d6ab1dc0SHong Zhang while (i--) { 1845c280ed6aSJed Brown PetscMPIInt icompleted; 1846c280ed6aSJed Brown ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 1847d6ab1dc0SHong Zhang } 1848d6ab1dc0SHong Zhang ierr = PetscFree(rwaits);CHKERRQ(ierr); 1849d6ab1dc0SHong Zhang if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);} 1850d6ab1dc0SHong Zhang ierr = PetscFree(len_si);CHKERRQ(ierr); 1851d6ab1dc0SHong Zhang ierr = PetscFree(len_ri);CHKERRQ(ierr); 1852d6ab1dc0SHong Zhang ierr = PetscFree(swaits);CHKERRQ(ierr); 1853d6ab1dc0SHong Zhang ierr = PetscFree(sstatus);CHKERRQ(ierr); 1854d6ab1dc0SHong Zhang ierr = PetscFree(buf_s);CHKERRQ(ierr); 1855d6ab1dc0SHong Zhang 1856d6ab1dc0SHong Zhang /* compute the local portion of C (mpi mat) */ 1857d6ab1dc0SHong Zhang /*------------------------------------------*/ 1858d6ab1dc0SHong Zhang /* allocate bi array and free space for accumulating nonzero column info */ 1859854ce69bSBarry Smith ierr = PetscMalloc1(pn+1,&bi);CHKERRQ(ierr); 1860d6ab1dc0SHong Zhang bi[0] = 0; 1861d6ab1dc0SHong Zhang 1862d6ab1dc0SHong Zhang /* set initial free space to be fill*(nnz(P) + nnz(AP)) */ 1863d6ab1dc0SHong Zhang nnz = fill*(pdti[pn] + poti[pon] + ai[am]); 1864a2ea699eSBarry Smith ierr = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr); 1865d6ab1dc0SHong Zhang current_space = free_space; 1866d6ab1dc0SHong Zhang 1867dcca6d9dSJed Brown ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr); 1868d6ab1dc0SHong Zhang for (k=0; k<merge->nrecv; k++) { 1869d6ab1dc0SHong Zhang buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 1870d6ab1dc0SHong Zhang nrows = *buf_ri_k[k]; 1871d6ab1dc0SHong Zhang nextrow[k] = buf_ri_k[k] + 1; /* next row number of k-th recved i-structure */ 18722205254eSKarl Rupp nextci[k] = buf_ri_k[k] + (nrows + 1); /* points to the next i-structure of k-th recieved i-structure */ 1873d6ab1dc0SHong Zhang } 1874d6ab1dc0SHong Zhang 1875d6ab1dc0SHong Zhang ierr = MatPreallocateInitialize(comm,pn,A->cmap->n,dnz,onz);CHKERRQ(ierr); 1876d6ab1dc0SHong Zhang rmax = 0; 1877d6ab1dc0SHong Zhang for (i=0; i<pn; i++) { 1878d6ab1dc0SHong Zhang /* add pdt[i,:]*AP into lnk */ 1879d6ab1dc0SHong Zhang pnz = pdti[i+1] - pdti[i]; 1880d6ab1dc0SHong Zhang ptJ = pdtj + pdti[i]; 1881d6ab1dc0SHong Zhang for (j=0; j<pnz; j++) { 1882d6ab1dc0SHong Zhang row = ptJ[j]; /* row of AP == col of Pt */ 1883d6ab1dc0SHong Zhang anz = ai[row+1] - ai[row]; 1884d6ab1dc0SHong Zhang Jptr = aj + ai[row]; 1885d6ab1dc0SHong Zhang /* add non-zero cols of AP into the sorted linked list lnk */ 1886d6ab1dc0SHong Zhang ierr = PetscLLCondensedAddSorted_Scalable(anz,Jptr,lnk);CHKERRQ(ierr); 1887d6ab1dc0SHong Zhang } 1888d6ab1dc0SHong Zhang 1889d6ab1dc0SHong Zhang /* add received col data into lnk */ 1890d6ab1dc0SHong Zhang for (k=0; k<merge->nrecv; k++) { /* k-th received message */ 1891d6ab1dc0SHong Zhang if (i == *nextrow[k]) { /* i-th row */ 1892d6ab1dc0SHong Zhang nzi = *(nextci[k]+1) - *nextci[k]; 1893d6ab1dc0SHong Zhang Jptr = buf_rj[k] + *nextci[k]; 1894d6ab1dc0SHong Zhang ierr = PetscLLCondensedAddSorted_Scalable(nzi,Jptr,lnk);CHKERRQ(ierr); 1895d6ab1dc0SHong Zhang nextrow[k]++; nextci[k]++; 1896d6ab1dc0SHong Zhang } 1897d6ab1dc0SHong Zhang } 1898d6ab1dc0SHong Zhang nnz = lnk[0]; 1899d6ab1dc0SHong Zhang 1900d6ab1dc0SHong Zhang /* if free space is not available, make more free space */ 1901d6ab1dc0SHong Zhang if (current_space->local_remaining<nnz) { 1902d6ab1dc0SHong Zhang ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 1903d6ab1dc0SHong Zhang nspacedouble++; 1904d6ab1dc0SHong Zhang } 1905d6ab1dc0SHong Zhang /* copy data into free space, then initialize lnk */ 1906d6ab1dc0SHong Zhang ierr = PetscLLCondensedClean_Scalable(nnz,current_space->array,lnk);CHKERRQ(ierr); 1907d6ab1dc0SHong Zhang ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr); 19082205254eSKarl Rupp 1909d6ab1dc0SHong Zhang current_space->array += nnz; 1910d6ab1dc0SHong Zhang current_space->local_used += nnz; 1911d6ab1dc0SHong Zhang current_space->local_remaining -= nnz; 19122205254eSKarl Rupp 1913d6ab1dc0SHong Zhang bi[i+1] = bi[i] + nnz; 1914d6ab1dc0SHong Zhang if (nnz > rmax) rmax = nnz; 1915d6ab1dc0SHong Zhang } 1916d6ab1dc0SHong Zhang ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr); 1917d6ab1dc0SHong Zhang 1918854ce69bSBarry Smith ierr = PetscMalloc1(bi[pn]+1,&bj);CHKERRQ(ierr); 1919d6ab1dc0SHong Zhang ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr); 1920118727c9SMark F. Adams afill_tmp = (PetscReal)bi[pn]/(pdti[pn] + poti[pon] + ai[am]+1); 1921d6ab1dc0SHong Zhang if (afill_tmp > afill) afill = afill_tmp; 1922d6ab1dc0SHong Zhang ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr); 1923d6ab1dc0SHong Zhang ierr = MatDestroy(&POt);CHKERRQ(ierr); 1924d6ab1dc0SHong Zhang ierr = MatDestroy(&PDt);CHKERRQ(ierr); 1925d6ab1dc0SHong Zhang 1926d6ab1dc0SHong Zhang /* create symbolic parallel matrix Cmpi - why cannot be assembled in Numeric part */ 1927d6ab1dc0SHong Zhang /*----------------------------------------------------------------------------------*/ 19281795a4d1SJed Brown ierr = PetscCalloc1(rmax+1,&vals);CHKERRQ(ierr); 1929d6ab1dc0SHong Zhang 1930d6ab1dc0SHong Zhang ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr); 1931d6ab1dc0SHong Zhang ierr = MatSetSizes(Cmpi,pn,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 193233d57670SJed Brown ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(A->cmap->bs));CHKERRQ(ierr); 1933d6ab1dc0SHong Zhang ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr); 1934d6ab1dc0SHong Zhang ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr); 1935d6ab1dc0SHong Zhang ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 1936d6ab1dc0SHong Zhang ierr = MatSetBlockSize(Cmpi,1);CHKERRQ(ierr); 1937d6ab1dc0SHong Zhang for (i=0; i<pn; i++) { 1938d6ab1dc0SHong Zhang row = i + rstart; 1939d6ab1dc0SHong Zhang nnz = bi[i+1] - bi[i]; 1940d6ab1dc0SHong Zhang Jptr = bj + bi[i]; 1941d6ab1dc0SHong Zhang ierr = MatSetValues(Cmpi,1,&row,nnz,Jptr,vals,INSERT_VALUES);CHKERRQ(ierr); 1942d6ab1dc0SHong Zhang } 1943d6ab1dc0SHong Zhang ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1944d6ab1dc0SHong Zhang ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1945d6ab1dc0SHong Zhang ierr = PetscFree(vals);CHKERRQ(ierr); 1946d6ab1dc0SHong Zhang 1947d6ab1dc0SHong Zhang merge->bi = bi; 1948d6ab1dc0SHong Zhang merge->bj = bj; 1949d6ab1dc0SHong Zhang merge->coi = coi; 1950d6ab1dc0SHong Zhang merge->coj = coj; 1951d6ab1dc0SHong Zhang merge->buf_ri = buf_ri; 1952d6ab1dc0SHong Zhang merge->buf_rj = buf_rj; 1953d6ab1dc0SHong Zhang merge->owners_co = owners_co; 1954d6ab1dc0SHong Zhang merge->destroy = Cmpi->ops->destroy; 1955d6ab1dc0SHong Zhang merge->duplicate = Cmpi->ops->duplicate; 1956d6ab1dc0SHong Zhang 19576da69ca6SHong Zhang Cmpi->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ; 1958d6ab1dc0SHong Zhang Cmpi->ops->destroy = MatDestroy_MPIAIJ_PtAP; 1959c9ba551fSBarry Smith Cmpi->ops->duplicate = MatDuplicate_MPIAIJ_MatPtAP; 1960d6ab1dc0SHong Zhang 1961d6ab1dc0SHong Zhang /* attach the supporting struct to Cmpi for reuse */ 1962d6ab1dc0SHong Zhang c = (Mat_MPIAIJ*)Cmpi->data; 19632205254eSKarl Rupp 1964d6ab1dc0SHong Zhang c->ptap = ptap; 19650298fd71SBarry Smith ptap->api = NULL; 19660298fd71SBarry Smith ptap->apj = NULL; 1967d6ab1dc0SHong Zhang ptap->merge = merge; 1968d6ab1dc0SHong Zhang ptap->rmax = rmax; 19690298fd71SBarry Smith ptap->apa = NULL; 1970f96d37fcSHong Zhang 1971f96d37fcSHong Zhang *C = Cmpi; 1972f96d37fcSHong Zhang #if defined(PETSC_USE_INFO) 1973f96d37fcSHong Zhang if (bi[pn] != 0) { 197457622a8eSBarry Smith ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr); 197557622a8eSBarry Smith ierr = PetscInfo1(Cmpi,"Use MatTransposeMatMult(A,B,MatReuse,%g,&C) for best performance.\n",(double)afill);CHKERRQ(ierr); 1976f96d37fcSHong Zhang } else { 1977f96d37fcSHong Zhang ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr); 1978f96d37fcSHong Zhang } 1979f96d37fcSHong Zhang #endif 1980f96d37fcSHong Zhang PetscFunctionReturn(0); 1981f96d37fcSHong Zhang } 1982