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 135e5acdf2Sstefano_zampini #if defined(PETSC_HAVE_HYPRE) 145e5acdf2Sstefano_zampini PETSC_INTERN PetscErrorCode MatMatMultSymbolic_AIJ_AIJ_wHYPRE(Mat,Mat,PetscReal,Mat*); 155e5acdf2Sstefano_zampini #endif 165e5acdf2Sstefano_zampini 17150d2497SBarry Smith PETSC_INTERN PetscErrorCode MatMatMult_MPIAIJ_MPIAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill, Mat *C) 182499ec78SHong Zhang { 192499ec78SHong Zhang PetscErrorCode ierr; 205e5acdf2Sstefano_zampini #if defined(PETSC_HAVE_HYPRE) 21a9d7e0ccSandi selinger const char *algTypes[4] = {"scalable","nonscalable","seqmpi","hypre"}; 22a9d7e0ccSandi selinger PetscInt nalg = 4; 235e5acdf2Sstefano_zampini #else 24a9d7e0ccSandi selinger const char *algTypes[3] = {"scalable","nonscalable","seqmpi"}; 25a9d7e0ccSandi selinger PetscInt nalg = 3; 265e5acdf2Sstefano_zampini #endif 2760106aa7SHong Zhang PetscInt alg = 1; /* set nonscalable algorithm as default */ 28e55be12dSHong Zhang MPI_Comm comm; 295e68f438SHong Zhang PetscBool flg; 302499ec78SHong Zhang 312499ec78SHong Zhang PetscFunctionBegin; 322499ec78SHong Zhang if (scall == MAT_INITIAL_MATRIX) { 33e55be12dSHong Zhang ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 346c4ed002SBarry Smith if (A->cmap->rstart != B->rmap->rstart || A->cmap->rend != B->rmap->rend) SETERRQ4(comm,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, (%D, %D) != (%D,%D)",A->cmap->rstart,A->cmap->rend,B->rmap->rstart,B->rmap->rend); 35e55be12dSHong Zhang 36715a5346SHong Zhang ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"MatMatMult","Mat");CHKERRQ(ierr); 375e68f438SHong Zhang ierr = PetscOptionsEList("-matmatmult_via","Algorithmic approach","MatMatMult",algTypes,nalg,algTypes[1],&alg,&flg);CHKERRQ(ierr); 380fc8cf34SHong Zhang ierr = PetscOptionsEnd();CHKERRQ(ierr); 390fc8cf34SHong Zhang 4060106aa7SHong Zhang if (!flg && B->cmap->N > 100000) { /* may switch to scalable algorithm as default */ 415e68f438SHong Zhang MatInfo Ainfo,Binfo; 4260106aa7SHong Zhang PetscInt nz_local; 4360106aa7SHong Zhang PetscBool alg_scalable_loc=PETSC_FALSE,alg_scalable; 445e68f438SHong Zhang 455e68f438SHong Zhang ierr = MatGetInfo(A,MAT_LOCAL,&Ainfo);CHKERRQ(ierr); 465e68f438SHong Zhang ierr = MatGetInfo(B,MAT_LOCAL,&Binfo);CHKERRQ(ierr); 475e68f438SHong Zhang nz_local = (PetscInt)(Ainfo.nz_allocated + Binfo.nz_allocated); 485e68f438SHong Zhang 4960106aa7SHong Zhang if (B->cmap->N > fill*nz_local) alg_scalable_loc = PETSC_TRUE; 5060106aa7SHong Zhang ierr = MPIU_Allreduce(&alg_scalable_loc,&alg_scalable,1,MPIU_BOOL,MPI_LOR,comm);CHKERRQ(ierr); 5160106aa7SHong Zhang 5260106aa7SHong Zhang if (alg_scalable) { 5360106aa7SHong Zhang alg = 0; /* scalable algorithm would 50% slower than nonscalable algorithm */ 5460106aa7SHong Zhang ierr = PetscInfo2(B,"Use scalable algorithm, BN %D, fill*nz_allocated %g\n",B->cmap->N,fill*nz_local);CHKERRQ(ierr); 5560106aa7SHong Zhang } 565e68f438SHong Zhang } 575e68f438SHong Zhang 583ff4c91cSHong Zhang ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr); 590fc8cf34SHong Zhang switch (alg) { 600fc8cf34SHong Zhang case 1: 610fc8cf34SHong Zhang ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(A,B,fill,C);CHKERRQ(ierr); 620fc8cf34SHong Zhang break; 635e5acdf2Sstefano_zampini case 2: 64a9d7e0ccSandi selinger ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ_seqMPI(A,B,fill,C);CHKERRQ(ierr); 65a9d7e0ccSandi selinger break; 66a9d7e0ccSandi selinger #if defined(PETSC_HAVE_HYPRE) 67a9d7e0ccSandi selinger case 3: 685e5acdf2Sstefano_zampini ierr = MatMatMultSymbolic_AIJ_AIJ_wHYPRE(A,B,fill,C);CHKERRQ(ierr); 695e5acdf2Sstefano_zampini break; 705e5acdf2Sstefano_zampini #endif 710fc8cf34SHong Zhang default: 72b2405163SHong Zhang ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ(A,B,fill,C);CHKERRQ(ierr); 730fc8cf34SHong Zhang break; 740fc8cf34SHong Zhang } 753ff4c91cSHong Zhang ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr); 76*7d0a89b7SHong Zhang 77*7d0a89b7SHong Zhang if (alg == 0 || alg == 1) { 78*7d0a89b7SHong Zhang Mat_MPIAIJ *c = (Mat_MPIAIJ*)(*C)->data; 79*7d0a89b7SHong Zhang Mat_APMPI *ap = c->ap; 80*7d0a89b7SHong Zhang ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)(*C)),((PetscObject)(*C))->prefix,"MatFreeIntermediateDataStructures","Mat");CHKERRQ(ierr); 81*7d0a89b7SHong Zhang ap->freestruct = PETSC_FALSE; 82*7d0a89b7SHong Zhang ierr = PetscOptionsBool("-mat_freeintermediatedatastructures","Free intermediate data structures", "MatFreeIntermediateDataStructures",ap->freestruct,&ap->freestruct, NULL);CHKERRQ(ierr); 83*7d0a89b7SHong Zhang ierr = PetscOptionsEnd();CHKERRQ(ierr); 842499ec78SHong Zhang } 85*7d0a89b7SHong Zhang } 86*7d0a89b7SHong Zhang 873ff4c91cSHong Zhang ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr); 88598bc09dSHong Zhang ierr = (*(*C)->ops->matmultnumeric)(A,B,*C);CHKERRQ(ierr); 893ff4c91cSHong Zhang ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr); 902499ec78SHong Zhang PetscFunctionReturn(0); 912499ec78SHong Zhang } 922499ec78SHong Zhang 93a1a4e84aSHong Zhang PetscErrorCode MatDestroy_MPIAIJ_MatMatMult(Mat A) 94598bc09dSHong Zhang { 95598bc09dSHong Zhang PetscErrorCode ierr; 96598bc09dSHong Zhang Mat_MPIAIJ *a = (Mat_MPIAIJ*)A->data; 973cdca5ebSHong Zhang Mat_APMPI *ptap = a->ap; 98598bc09dSHong Zhang 99598bc09dSHong Zhang PetscFunctionBegin; 100b7f45c76SHong Zhang ierr = PetscFree2(ptap->startsj_s,ptap->startsj_r);CHKERRQ(ierr); 101598bc09dSHong Zhang ierr = PetscFree(ptap->bufa);CHKERRQ(ierr); 102a1a4e84aSHong Zhang ierr = MatDestroy(&ptap->P_loc);CHKERRQ(ierr); 103a1a4e84aSHong Zhang ierr = MatDestroy(&ptap->P_oth);CHKERRQ(ierr); 104ab1b013aSHong Zhang ierr = MatDestroy(&ptap->Pt);CHKERRQ(ierr); 105a1a4e84aSHong Zhang ierr = PetscFree(ptap->api);CHKERRQ(ierr); 106a1a4e84aSHong Zhang ierr = PetscFree(ptap->apj);CHKERRQ(ierr); 107598bc09dSHong Zhang ierr = PetscFree(ptap->apa);CHKERRQ(ierr); 108598bc09dSHong Zhang ierr = ptap->destroy(A);CHKERRQ(ierr); 109598bc09dSHong Zhang ierr = PetscFree(ptap);CHKERRQ(ierr); 110598bc09dSHong Zhang PetscFunctionReturn(0); 111598bc09dSHong Zhang } 112598bc09dSHong Zhang 113f8efd71cSHong Zhang PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable(Mat A,Mat P,Mat C) 114598bc09dSHong Zhang { 115598bc09dSHong Zhang PetscErrorCode ierr; 1164ae0847bSHong Zhang Mat_MPIAIJ *a =(Mat_MPIAIJ*)A->data,*c=(Mat_MPIAIJ*)C->data; 117598bc09dSHong Zhang Mat_SeqAIJ *ad =(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data; 118598bc09dSHong Zhang Mat_SeqAIJ *cd =(Mat_SeqAIJ*)(c->A)->data,*co=(Mat_SeqAIJ*)(c->B)->data; 119c12557a1SHong Zhang PetscScalar *cda=cd->a,*coa=co->a; 120598bc09dSHong Zhang Mat_SeqAIJ *p_loc,*p_oth; 1219ce11a7cSHong Zhang PetscScalar *apa,*ca; 122c12557a1SHong Zhang PetscInt cm =C->rmap->n; 1233cdca5ebSHong Zhang Mat_APMPI *ptap=c->ap; 124c12557a1SHong Zhang PetscInt *api,*apj,*apJ,i,k; 12529825010SHong Zhang PetscInt cstart=C->cmap->rstart; 126598bc09dSHong Zhang PetscInt cdnz,conz,k0,k1; 1279467f45fSHong Zhang MPI_Comm comm; 1289467f45fSHong Zhang PetscMPIInt size; 129598bc09dSHong Zhang 130598bc09dSHong Zhang PetscFunctionBegin; 1319467f45fSHong Zhang ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 1329467f45fSHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1339467f45fSHong Zhang 134*7d0a89b7SHong Zhang if (!ptap) SETERRQ(comm,PETSC_ERR_ARG_WRONGSTATE,"AP cannot be reused. Do not call MatFreeIntermediateDataStructures() or use '-mat_freeintermediatedatastructures'"); 13580da8df7SHong Zhang 136a1a4e84aSHong Zhang /* 1) get P_oth = ptap->P_oth and P_loc = ptap->P_loc */ 137598bc09dSHong Zhang /*-----------------------------------------------------*/ 138a1a4e84aSHong Zhang /* update numerical values of P_oth and P_loc */ 139b7f45c76SHong Zhang ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 140a1a4e84aSHong Zhang ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 141598bc09dSHong Zhang 142598bc09dSHong Zhang /* 2) compute numeric C_loc = A_loc*P = Ad*P_loc + Ao*P_oth */ 143598bc09dSHong Zhang /*----------------------------------------------------------*/ 144598bc09dSHong Zhang /* get data from symbolic products */ 145a1a4e84aSHong Zhang p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 146c12557a1SHong Zhang p_oth = NULL; 1479467f45fSHong Zhang if (size >1) { 1489467f45fSHong Zhang p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 1499467f45fSHong Zhang } 150598bc09dSHong Zhang 151598bc09dSHong Zhang /* get apa for storing dense row A[i,:]*P */ 152598bc09dSHong Zhang apa = ptap->apa; 153598bc09dSHong Zhang 15429825010SHong Zhang api = ptap->api; 15529825010SHong Zhang apj = ptap->apj; 156598bc09dSHong Zhang for (i=0; i<cm; i++) { 157c12557a1SHong Zhang /* compute apa = A[i,:]*P */ 158e497d3c8SHong Zhang AProw_nonscalable(i,ad,ao,p_loc,p_oth,apa); 159598bc09dSHong Zhang 160598bc09dSHong Zhang /* set values in C */ 161598bc09dSHong Zhang apJ = apj + api[i]; 162598bc09dSHong Zhang cdnz = cd->i[i+1] - cd->i[i]; 163598bc09dSHong Zhang conz = co->i[i+1] - co->i[i]; 164598bc09dSHong Zhang 165598bc09dSHong Zhang /* 1st off-diagoanl part of C */ 166598bc09dSHong Zhang ca = coa + co->i[i]; 167598bc09dSHong Zhang k = 0; 168598bc09dSHong Zhang for (k0=0; k0<conz; k0++) { 169598bc09dSHong Zhang if (apJ[k] >= cstart) break; 170598bc09dSHong Zhang ca[k0] = apa[apJ[k]]; 1715cab4f04SHong Zhang apa[apJ[k++]] = 0.0; 172598bc09dSHong Zhang } 173598bc09dSHong Zhang 174598bc09dSHong Zhang /* diagonal part of C */ 175598bc09dSHong Zhang ca = cda + cd->i[i]; 176598bc09dSHong Zhang for (k1=0; k1<cdnz; k1++) { 177598bc09dSHong Zhang ca[k1] = apa[apJ[k]]; 1785cab4f04SHong Zhang apa[apJ[k++]] = 0.0; 179598bc09dSHong Zhang } 180598bc09dSHong Zhang 181598bc09dSHong Zhang /* 2nd off-diagoanl part of C */ 182598bc09dSHong Zhang ca = coa + co->i[i]; 183598bc09dSHong Zhang for (; k0<conz; k0++) { 184598bc09dSHong Zhang ca[k0] = apa[apJ[k]]; 1855cab4f04SHong Zhang apa[apJ[k++]] = 0.0; 186598bc09dSHong Zhang } 187598bc09dSHong Zhang } 188598bc09dSHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 189598bc09dSHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 19080da8df7SHong Zhang 191*7d0a89b7SHong Zhang if (ptap->freestruct) { 19280da8df7SHong Zhang ierr = MatFreeIntermediateDataStructures(C);CHKERRQ(ierr); 19380da8df7SHong Zhang } 194598bc09dSHong Zhang PetscFunctionReturn(0); 195598bc09dSHong Zhang } 196598bc09dSHong Zhang 1970fc8cf34SHong Zhang PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat A,Mat P,PetscReal fill,Mat *C) 198598bc09dSHong Zhang { 199598bc09dSHong Zhang PetscErrorCode ierr; 200ce94432eSBarry Smith MPI_Comm comm; 2019467f45fSHong Zhang PetscMPIInt size; 202bfcd1a73SHong Zhang Mat Cmpi; 2033cdca5ebSHong Zhang Mat_APMPI *ptap; 2040298fd71SBarry Smith PetscFreeSpaceList free_space=NULL,current_space=NULL; 2054ae0847bSHong Zhang Mat_MPIAIJ *a =(Mat_MPIAIJ*)A->data,*c; 206bfcd1a73SHong Zhang Mat_SeqAIJ *ad =(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_loc,*p_oth; 2077dbaa2c6Sandi selinger PetscInt *pi_loc,*pj_loc,*pi_oth,*pj_oth,*dnz,*onz; 2084ae0847bSHong Zhang PetscInt *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,rstart=A->rmap->rstart; 209d14fa243SHong Zhang PetscInt *lnk,i,pnz,row,*api,*apj,*Jptr,apnz,nspacedouble=0,j,nzi; 2105cab4f04SHong Zhang PetscInt am=A->rmap->n,pN=P->cmap->N,pn=P->cmap->n,pm=P->rmap->n; 211598bc09dSHong Zhang PetscBT lnkbt; 212598bc09dSHong Zhang PetscScalar *apa; 213bfcd1a73SHong Zhang PetscReal afill; 214b92f168fSBarry Smith MatType mtype; 215598bc09dSHong Zhang 216598bc09dSHong Zhang PetscFunctionBegin; 217ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 2189467f45fSHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 2199467f45fSHong Zhang 2203cdca5ebSHong Zhang /* create struct Mat_APMPI and attached it to C later */ 221b00a9115SJed Brown ierr = PetscNew(&ptap);CHKERRQ(ierr); 222598bc09dSHong Zhang 223598bc09dSHong Zhang /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */ 224b7f45c76SHong Zhang ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 22519f0e57aSHong Zhang 226598bc09dSHong Zhang /* get P_loc by taking all local rows of P */ 227a1a4e84aSHong Zhang ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 228598bc09dSHong Zhang 229a1a4e84aSHong Zhang p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 230598bc09dSHong Zhang pi_loc = p_loc->i; pj_loc = p_loc->j; 2319467f45fSHong Zhang if (size > 1) { 2329467f45fSHong Zhang p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 233598bc09dSHong Zhang pi_oth = p_oth->i; pj_oth = p_oth->j; 23420e3dc75SHong Zhang } else { 23520e3dc75SHong Zhang p_oth = NULL; 23620e3dc75SHong Zhang pi_oth = NULL; pj_oth = NULL; 2379467f45fSHong Zhang } 238598bc09dSHong Zhang 239598bc09dSHong Zhang /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */ 240598bc09dSHong Zhang /*-------------------------------------------------------------------*/ 241854ce69bSBarry Smith ierr = PetscMalloc1(am+2,&api);CHKERRQ(ierr); 242a1a4e84aSHong Zhang ptap->api = api; 243598bc09dSHong Zhang api[0] = 0; 244598bc09dSHong Zhang 2455cab4f04SHong Zhang /* create and initialize a linked list */ 2465cab4f04SHong Zhang ierr = PetscLLCondensedCreate(pN,pN,&lnk,&lnkbt);CHKERRQ(ierr); 247598bc09dSHong Zhang 248bfcd1a73SHong Zhang /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */ 249f91af8c7SBarry Smith ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(adi[am],PetscIntSumTruncate(aoi[am],pi_loc[pm]))),&free_space);CHKERRQ(ierr); 250598bc09dSHong Zhang current_space = free_space; 251598bc09dSHong Zhang 252598bc09dSHong Zhang ierr = MatPreallocateInitialize(comm,am,pn,dnz,onz);CHKERRQ(ierr); 253598bc09dSHong Zhang for (i=0; i<am; i++) { 254598bc09dSHong Zhang /* diagonal portion of A */ 255598bc09dSHong Zhang nzi = adi[i+1] - adi[i]; 256598bc09dSHong Zhang for (j=0; j<nzi; j++) { 257598bc09dSHong Zhang row = *adj++; 258598bc09dSHong Zhang pnz = pi_loc[row+1] - pi_loc[row]; 259598bc09dSHong Zhang Jptr = pj_loc + pi_loc[row]; 260598bc09dSHong Zhang /* add non-zero cols of P into the sorted linked list lnk */ 2611fbbb359SHong Zhang ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 262598bc09dSHong Zhang } 263598bc09dSHong Zhang /* off-diagonal portion of A */ 264598bc09dSHong Zhang nzi = aoi[i+1] - aoi[i]; 265598bc09dSHong Zhang for (j=0; j<nzi; j++) { 266598bc09dSHong Zhang row = *aoj++; 267598bc09dSHong Zhang pnz = pi_oth[row+1] - pi_oth[row]; 268598bc09dSHong Zhang Jptr = pj_oth + pi_oth[row]; 2691fbbb359SHong Zhang ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 270598bc09dSHong Zhang } 271598bc09dSHong Zhang 272d14fa243SHong Zhang apnz = lnk[0]; 273598bc09dSHong Zhang api[i+1] = api[i] + apnz; 274598bc09dSHong Zhang 275598bc09dSHong Zhang /* if free space is not available, double the total space in the list */ 276598bc09dSHong Zhang if (current_space->local_remaining<apnz) { 277f91af8c7SBarry Smith ierr = PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),¤t_space);CHKERRQ(ierr); 278598bc09dSHong Zhang nspacedouble++; 279598bc09dSHong Zhang } 280598bc09dSHong Zhang 281598bc09dSHong Zhang /* Copy data into free space, then initialize lnk */ 2821fbbb359SHong Zhang ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 283598bc09dSHong Zhang ierr = MatPreallocateSet(i+rstart,apnz,current_space->array,dnz,onz);CHKERRQ(ierr); 2842205254eSKarl Rupp 285598bc09dSHong Zhang current_space->array += apnz; 286598bc09dSHong Zhang current_space->local_used += apnz; 287598bc09dSHong Zhang current_space->local_remaining -= apnz; 288598bc09dSHong Zhang } 289598bc09dSHong Zhang 290598bc09dSHong Zhang /* Allocate space for apj, initialize apj, and */ 291598bc09dSHong Zhang /* destroy list of free space and other temporary array(s) */ 292854ce69bSBarry Smith ierr = PetscMalloc1(api[am]+1,&ptap->apj);CHKERRQ(ierr); 293a1a4e84aSHong Zhang apj = ptap->apj; 294a1a4e84aSHong Zhang ierr = PetscFreeSpaceContiguous(&free_space,ptap->apj);CHKERRQ(ierr); 295598bc09dSHong Zhang ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 296598bc09dSHong Zhang 297f84c536bSHong Zhang /* malloc apa to store dense row A[i,:]*P */ 2981795a4d1SJed Brown ierr = PetscCalloc1(pN,&apa);CHKERRQ(ierr); 2992205254eSKarl Rupp 300f84c536bSHong Zhang ptap->apa = apa; 301f84c536bSHong Zhang 302bfcd1a73SHong Zhang /* create and assemble symbolic parallel matrix Cmpi */ 303598bc09dSHong Zhang /*----------------------------------------------------*/ 304bfcd1a73SHong Zhang ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr); 305bfcd1a73SHong Zhang ierr = MatSetSizes(Cmpi,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 30633d57670SJed Brown ierr = MatSetBlockSizesFromMats(Cmpi,A,P);CHKERRQ(ierr); 307a2f3521dSMark F. Adams 308b92f168fSBarry Smith ierr = MatGetType(A,&mtype);CHKERRQ(ierr); 309b92f168fSBarry Smith ierr = MatSetType(Cmpi,mtype);CHKERRQ(ierr); 310bfcd1a73SHong Zhang ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr); 31117f75172Sandi selinger 312904d1e70Sandi selinger ierr = MatSetValues_MPIAIJ_CopyFromCSRFormat_Symbolic(Cmpi, apj, api);CHKERRQ(ierr); 313bfcd1a73SHong Zhang ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 314bfcd1a73SHong Zhang ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 3157dbaa2c6Sandi selinger ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 316598bc09dSHong Zhang 317bfcd1a73SHong Zhang ptap->destroy = Cmpi->ops->destroy; 318bfcd1a73SHong Zhang ptap->duplicate = Cmpi->ops->duplicate; 319f8efd71cSHong Zhang Cmpi->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable; 320bfcd1a73SHong Zhang Cmpi->ops->destroy = MatDestroy_MPIAIJ_MatMatMult; 3214624976aSHong Zhang Cmpi->ops->freeintermediatedatastructures = MatFreeIntermediateDataStructures_MPIAIJ_AP; 322598bc09dSHong Zhang 323bfcd1a73SHong Zhang /* attach the supporting struct to Cmpi for reuse */ 324bfcd1a73SHong Zhang c = (Mat_MPIAIJ*)Cmpi->data; 3253cdca5ebSHong Zhang c->ap = ptap; 326598bc09dSHong Zhang 327bfcd1a73SHong Zhang *C = Cmpi; 328bfcd1a73SHong Zhang 329bfcd1a73SHong Zhang /* set MatInfo */ 330118727c9SMark F. Adams afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1) + 1.e-5; 331bfcd1a73SHong Zhang if (afill < 1.0) afill = 1.0; 332bfcd1a73SHong Zhang Cmpi->info.mallocs = nspacedouble; 333bfcd1a73SHong Zhang Cmpi->info.fill_ratio_given = fill; 334bfcd1a73SHong Zhang Cmpi->info.fill_ratio_needed = afill; 335bfcd1a73SHong Zhang 336bfcd1a73SHong Zhang #if defined(PETSC_USE_INFO) 337bfcd1a73SHong Zhang if (api[am]) { 33857622a8eSBarry Smith ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr); 33957622a8eSBarry Smith ierr = PetscInfo1(Cmpi,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr); 340bfcd1a73SHong Zhang } else { 341bfcd1a73SHong Zhang ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr); 342bfcd1a73SHong Zhang } 343bfcd1a73SHong Zhang #endif 344598bc09dSHong Zhang PetscFunctionReturn(0); 345598bc09dSHong Zhang } 346598bc09dSHong Zhang 347150d2497SBarry Smith PETSC_INTERN PetscErrorCode MatMatMult_MPIAIJ_MPIDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 3489123193aSHong Zhang { 3499123193aSHong Zhang PetscErrorCode ierr; 3509123193aSHong Zhang 3519123193aSHong Zhang PetscFunctionBegin; 3529123193aSHong Zhang if (scall == MAT_INITIAL_MATRIX) { 3533ff4c91cSHong Zhang ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr); 3549123193aSHong Zhang ierr = MatMatMultSymbolic_MPIAIJ_MPIDense(A,B,fill,C);CHKERRQ(ierr); 3553ff4c91cSHong Zhang ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr); 3569123193aSHong Zhang } 3573ff4c91cSHong Zhang ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr); 3589123193aSHong Zhang ierr = MatMatMultNumeric_MPIAIJ_MPIDense(A,B,*C);CHKERRQ(ierr); 3593ff4c91cSHong Zhang ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr); 3609123193aSHong Zhang PetscFunctionReturn(0); 3619123193aSHong Zhang } 3629123193aSHong Zhang 3634324174eSBarry Smith typedef struct { 3644324174eSBarry Smith Mat workB; 3654324174eSBarry Smith PetscScalar *rvalues,*svalues; 3664324174eSBarry Smith MPI_Request *rwaits,*swaits; 3674324174eSBarry Smith } MPIAIJ_MPIDense; 3684324174eSBarry Smith 36996b95a6bSBarry Smith PetscErrorCode MatMPIAIJ_MPIDenseDestroy(void *ctx) 3704324174eSBarry Smith { 3714324174eSBarry Smith MPIAIJ_MPIDense *contents = (MPIAIJ_MPIDense*) ctx; 3724324174eSBarry Smith PetscErrorCode ierr; 3734324174eSBarry Smith 3744324174eSBarry Smith PetscFunctionBegin; 3756bf464f9SBarry Smith ierr = MatDestroy(&contents->workB);CHKERRQ(ierr); 3764324174eSBarry Smith ierr = PetscFree4(contents->rvalues,contents->svalues,contents->rwaits,contents->swaits);CHKERRQ(ierr); 3774324174eSBarry Smith ierr = PetscFree(contents);CHKERRQ(ierr); 3784324174eSBarry Smith PetscFunctionReturn(0); 3794324174eSBarry Smith } 3804324174eSBarry Smith 3816eb45d04SBarry Smith #include <petsc/private/vecscatterimpl.h> 382fd4e9aacSBarry Smith /* 383fd4e9aacSBarry Smith This is a "dummy function" that handles the case where matrix C was created as a dense matrix 384fd4e9aacSBarry Smith directly by the user and passed to MatMatMult() with the MAT_REUSE_MATRIX option 385fd4e9aacSBarry Smith 386fd4e9aacSBarry Smith It is the same as MatMatMultSymbolic_MPIAIJ_MPIDense() except does not create C 3876eb45d04SBarry Smith 3886eb45d04SBarry Smith Developer Notes: This directly accesses information inside the VecScatter associated with the matrix-vector product 3896eb45d04SBarry Smith for this matrix. This is not desirable.. 3906eb45d04SBarry Smith 391fd4e9aacSBarry Smith */ 392fd4e9aacSBarry Smith PetscErrorCode MatMatMultNumeric_MPIDense(Mat A,Mat B,Mat C) 393fd4e9aacSBarry Smith { 394fd4e9aacSBarry Smith PetscErrorCode ierr; 395fd4e9aacSBarry Smith PetscBool flg; 396fd4e9aacSBarry Smith Mat_MPIAIJ *aij = (Mat_MPIAIJ*) A->data; 397fd4e9aacSBarry Smith PetscInt nz = aij->B->cmap->n; 398fd4e9aacSBarry Smith PetscContainer container; 399fd4e9aacSBarry Smith MPIAIJ_MPIDense *contents; 400fd4e9aacSBarry Smith VecScatter ctx = aij->Mvctx; 401fd4e9aacSBarry Smith VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata; 402fd4e9aacSBarry Smith VecScatter_MPI_General *to = (VecScatter_MPI_General*) ctx->todata; 403fd4e9aacSBarry Smith 404fd4e9aacSBarry Smith PetscFunctionBegin; 405fd4e9aacSBarry Smith ierr = PetscObjectTypeCompare((PetscObject)B,MATMPIDENSE,&flg);CHKERRQ(ierr); 406fd4e9aacSBarry Smith if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Second matrix must be mpidense"); 407fd4e9aacSBarry Smith 408fd4e9aacSBarry Smith /* Handle case where where user provided the final C matrix rather than calling MatMatMult() with MAT_INITIAL_MATRIX*/ 409fd4e9aacSBarry Smith ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&flg);CHKERRQ(ierr); 410fd4e9aacSBarry Smith if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"First matrix must be MPIAIJ"); 411fd4e9aacSBarry Smith 412fd4e9aacSBarry Smith C->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIDense; 413fd4e9aacSBarry Smith 414fd4e9aacSBarry Smith ierr = PetscNew(&contents);CHKERRQ(ierr); 415fd4e9aacSBarry Smith /* Create work matrix used to store off processor rows of B needed for local product */ 416fd4e9aacSBarry Smith ierr = MatCreateSeqDense(PETSC_COMM_SELF,nz,B->cmap->N,NULL,&contents->workB);CHKERRQ(ierr); 417fd4e9aacSBarry Smith /* Create work arrays needed */ 418fd4e9aacSBarry Smith ierr = PetscMalloc4(B->cmap->N*from->starts[from->n],&contents->rvalues, 419fd4e9aacSBarry Smith B->cmap->N*to->starts[to->n],&contents->svalues, 420fd4e9aacSBarry Smith from->n,&contents->rwaits, 421fd4e9aacSBarry Smith to->n,&contents->swaits);CHKERRQ(ierr); 422fd4e9aacSBarry Smith 423fd4e9aacSBarry Smith ierr = PetscContainerCreate(PetscObjectComm((PetscObject)A),&container);CHKERRQ(ierr); 424fd4e9aacSBarry Smith ierr = PetscContainerSetPointer(container,contents);CHKERRQ(ierr); 425fd4e9aacSBarry Smith ierr = PetscContainerSetUserDestroy(container,MatMPIAIJ_MPIDenseDestroy);CHKERRQ(ierr); 426fd4e9aacSBarry Smith ierr = PetscObjectCompose((PetscObject)C,"workB",(PetscObject)container);CHKERRQ(ierr); 427fd4e9aacSBarry Smith ierr = PetscContainerDestroy(&container);CHKERRQ(ierr); 428fd4e9aacSBarry Smith 429fd4e9aacSBarry Smith ierr = (*C->ops->matmultnumeric)(A,B,C);CHKERRQ(ierr); 430fd4e9aacSBarry Smith PetscFunctionReturn(0); 431fd4e9aacSBarry Smith } 432fd4e9aacSBarry Smith 4339123193aSHong Zhang PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIDense(Mat A,Mat B,PetscReal fill,Mat *C) 4349123193aSHong Zhang { 4359123193aSHong Zhang PetscErrorCode ierr; 436f32d5d43SBarry Smith Mat_MPIAIJ *aij = (Mat_MPIAIJ*) A->data; 437d0f46423SBarry Smith PetscInt nz = aij->B->cmap->n; 438bf0cc555SLisandro Dalcin PetscContainer container; 4394324174eSBarry Smith MPIAIJ_MPIDense *contents; 4404324174eSBarry Smith VecScatter ctx = aij->Mvctx; 4414324174eSBarry Smith VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata; 4424324174eSBarry Smith VecScatter_MPI_General *to = (VecScatter_MPI_General*) ctx->todata; 443d0f46423SBarry Smith PetscInt m = A->rmap->n,n=B->cmap->n; 4449123193aSHong Zhang 4459123193aSHong Zhang PetscFunctionBegin; 446ce94432eSBarry Smith ierr = MatCreate(PetscObjectComm((PetscObject)B),C);CHKERRQ(ierr); 447d0f46423SBarry Smith ierr = MatSetSizes(*C,m,n,A->rmap->N,B->cmap->N);CHKERRQ(ierr); 44833d57670SJed Brown ierr = MatSetBlockSizesFromMats(*C,A,B);CHKERRQ(ierr); 449cb2480eaSBarry Smith ierr = MatSetType(*C,MATMPIDENSE);CHKERRQ(ierr); 4500298fd71SBarry Smith ierr = MatMPIDenseSetPreallocation(*C,NULL);CHKERRQ(ierr); 451cb2480eaSBarry Smith ierr = MatAssemblyBegin(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 452cb2480eaSBarry Smith ierr = MatAssemblyEnd(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 4532205254eSKarl Rupp 454d73949e8SHong Zhang (*C)->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIDense; 455f32d5d43SBarry Smith 456b00a9115SJed Brown ierr = PetscNew(&contents);CHKERRQ(ierr); 457f32d5d43SBarry Smith /* Create work matrix used to store off processor rows of B needed for local product */ 4580298fd71SBarry Smith ierr = MatCreateSeqDense(PETSC_COMM_SELF,nz,B->cmap->N,NULL,&contents->workB);CHKERRQ(ierr); 4594324174eSBarry Smith /* Create work arrays needed */ 460dcca6d9dSJed Brown ierr = PetscMalloc4(B->cmap->N*from->starts[from->n],&contents->rvalues, 461dcca6d9dSJed Brown B->cmap->N*to->starts[to->n],&contents->svalues, 462dcca6d9dSJed Brown from->n,&contents->rwaits, 463dcca6d9dSJed Brown to->n,&contents->swaits);CHKERRQ(ierr); 4644324174eSBarry Smith 465ce94432eSBarry Smith ierr = PetscContainerCreate(PetscObjectComm((PetscObject)A),&container);CHKERRQ(ierr); 466bf0cc555SLisandro Dalcin ierr = PetscContainerSetPointer(container,contents);CHKERRQ(ierr); 46796b95a6bSBarry Smith ierr = PetscContainerSetUserDestroy(container,MatMPIAIJ_MPIDenseDestroy);CHKERRQ(ierr); 468bf0cc555SLisandro Dalcin ierr = PetscObjectCompose((PetscObject)(*C),"workB",(PetscObject)container);CHKERRQ(ierr); 469bf0cc555SLisandro Dalcin ierr = PetscContainerDestroy(&container);CHKERRQ(ierr); 4709123193aSHong Zhang PetscFunctionReturn(0); 4719123193aSHong Zhang } 4729123193aSHong Zhang 473f32d5d43SBarry Smith /* 4742636ff26SBarry Smith Performs an efficient scatter on the rows of B needed by this process; this is 4752636ff26SBarry Smith a modification of the VecScatterBegin_() routines. 476f32d5d43SBarry Smith */ 4774324174eSBarry Smith PetscErrorCode MatMPIDenseScatter(Mat A,Mat B,Mat C,Mat *outworkB) 478f32d5d43SBarry Smith { 479f32d5d43SBarry Smith Mat_MPIAIJ *aij = (Mat_MPIAIJ*)A->data; 480f32d5d43SBarry Smith PetscErrorCode ierr; 481f32d5d43SBarry Smith PetscScalar *b,*w,*svalues,*rvalues; 482f32d5d43SBarry Smith VecScatter ctx = aij->Mvctx; 483f32d5d43SBarry Smith VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata; 484f32d5d43SBarry Smith VecScatter_MPI_General *to = (VecScatter_MPI_General*) ctx->todata; 485f32d5d43SBarry Smith PetscInt i,j,k; 486aa5bb8c0SSatish Balay PetscInt *sindices,*sstarts,*rindices,*rstarts; 487aa5bb8c0SSatish Balay PetscMPIInt *sprocs,*rprocs,nrecvs; 488f32d5d43SBarry Smith MPI_Request *swaits,*rwaits; 489ce94432eSBarry Smith MPI_Comm comm; 490d0f46423SBarry Smith PetscMPIInt tag = ((PetscObject)ctx)->tag,ncols = B->cmap->N, nrows = aij->B->cmap->n,imdex,nrowsB = B->rmap->n; 491f32d5d43SBarry Smith MPI_Status status; 4924324174eSBarry Smith MPIAIJ_MPIDense *contents; 493bf0cc555SLisandro Dalcin PetscContainer container; 4944324174eSBarry Smith Mat workB; 495f32d5d43SBarry Smith 496f32d5d43SBarry Smith PetscFunctionBegin; 497ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 498bf0cc555SLisandro Dalcin ierr = PetscObjectQuery((PetscObject)C,"workB",(PetscObject*)&container);CHKERRQ(ierr); 49929825010SHong Zhang if (!container) SETERRQ(comm,PETSC_ERR_PLIB,"Container does not exist"); 500bf0cc555SLisandro Dalcin ierr = PetscContainerGetPointer(container,(void**)&contents);CHKERRQ(ierr); 5014324174eSBarry Smith 5024324174eSBarry Smith workB = *outworkB = contents->workB; 503cf1a0672SHong 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); 504f32d5d43SBarry Smith sindices = to->indices; 505f32d5d43SBarry Smith sstarts = to->starts; 506f32d5d43SBarry Smith sprocs = to->procs; 5074324174eSBarry Smith swaits = contents->swaits; 5084324174eSBarry Smith svalues = contents->svalues; 509f32d5d43SBarry Smith 510f32d5d43SBarry Smith rindices = from->indices; 511f32d5d43SBarry Smith rstarts = from->starts; 512f32d5d43SBarry Smith rprocs = from->procs; 5134324174eSBarry Smith rwaits = contents->rwaits; 5144324174eSBarry Smith rvalues = contents->rvalues; 515f32d5d43SBarry Smith 5168c778c55SBarry Smith ierr = MatDenseGetArray(B,&b);CHKERRQ(ierr); 5178c778c55SBarry Smith ierr = MatDenseGetArray(workB,&w);CHKERRQ(ierr); 518f32d5d43SBarry Smith 519f32d5d43SBarry Smith for (i=0; i<from->n; i++) { 520f32d5d43SBarry Smith ierr = MPI_Irecv(rvalues+ncols*rstarts[i],ncols*(rstarts[i+1]-rstarts[i]),MPIU_SCALAR,rprocs[i],tag,comm,rwaits+i);CHKERRQ(ierr); 521f32d5d43SBarry Smith } 522f32d5d43SBarry Smith 523f32d5d43SBarry Smith for (i=0; i<to->n; i++) { 524f32d5d43SBarry Smith /* pack a message at a time */ 525f32d5d43SBarry Smith for (j=0; j<sstarts[i+1]-sstarts[i]; j++) { 526f32d5d43SBarry Smith for (k=0; k<ncols; k++) { 527f32d5d43SBarry Smith svalues[ncols*(sstarts[i] + j) + k] = b[sindices[sstarts[i]+j] + nrowsB*k]; 5282636ff26SBarry Smith } 5292636ff26SBarry Smith } 530f32d5d43SBarry Smith ierr = MPI_Isend(svalues+ncols*sstarts[i],ncols*(sstarts[i+1]-sstarts[i]),MPIU_SCALAR,sprocs[i],tag,comm,swaits+i);CHKERRQ(ierr); 531f32d5d43SBarry Smith } 5322636ff26SBarry Smith 533f32d5d43SBarry Smith nrecvs = from->n; 534f32d5d43SBarry Smith while (nrecvs) { 535f32d5d43SBarry Smith ierr = MPI_Waitany(from->n,rwaits,&imdex,&status);CHKERRQ(ierr); 536f32d5d43SBarry Smith nrecvs--; 537f32d5d43SBarry Smith /* unpack a message at a time */ 538f32d5d43SBarry Smith for (j=0; j<rstarts[imdex+1]-rstarts[imdex]; j++) { 539f32d5d43SBarry Smith for (k=0; k<ncols; k++) { 540f32d5d43SBarry Smith w[rindices[rstarts[imdex]+j] + nrows*k] = rvalues[ncols*(rstarts[imdex] + j) + k]; 5412636ff26SBarry Smith } 5422636ff26SBarry Smith } 543f32d5d43SBarry Smith } 544cb9801acSJed Brown if (to->n) {ierr = MPI_Waitall(to->n,swaits,to->sstatus);CHKERRQ(ierr);} 545f32d5d43SBarry Smith 5468c778c55SBarry Smith ierr = MatDenseRestoreArray(B,&b);CHKERRQ(ierr); 5478c778c55SBarry Smith ierr = MatDenseRestoreArray(workB,&w);CHKERRQ(ierr); 548f32d5d43SBarry Smith ierr = MatAssemblyBegin(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 549f32d5d43SBarry Smith ierr = MatAssemblyEnd(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 550f32d5d43SBarry Smith PetscFunctionReturn(0); 551f32d5d43SBarry Smith } 552f32d5d43SBarry Smith extern PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat,Mat,Mat); 553f32d5d43SBarry Smith 5549123193aSHong Zhang PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIDense(Mat A,Mat B,Mat C) 5559123193aSHong Zhang { 5569123193aSHong Zhang PetscErrorCode ierr; 557f32d5d43SBarry Smith Mat_MPIAIJ *aij = (Mat_MPIAIJ*)A->data; 558f32d5d43SBarry Smith Mat_MPIDense *bdense = (Mat_MPIDense*)B->data; 559f32d5d43SBarry Smith Mat_MPIDense *cdense = (Mat_MPIDense*)C->data; 560f32d5d43SBarry Smith Mat workB; 5619123193aSHong Zhang 5629123193aSHong Zhang PetscFunctionBegin; 563f32d5d43SBarry Smith /* diagonal block of A times all local rows of B*/ 564f32d5d43SBarry Smith ierr = MatMatMultNumeric_SeqAIJ_SeqDense(aij->A,bdense->A,cdense->A);CHKERRQ(ierr); 565f32d5d43SBarry Smith 566f32d5d43SBarry Smith /* get off processor parts of B needed to complete the product */ 5674324174eSBarry Smith ierr = MatMPIDenseScatter(A,B,C,&workB);CHKERRQ(ierr); 568f32d5d43SBarry Smith 569f32d5d43SBarry Smith /* off-diagonal block of A times nonlocal rows of B */ 570f32d5d43SBarry Smith ierr = MatMatMultNumericAdd_SeqAIJ_SeqDense(aij->B,workB,cdense->A);CHKERRQ(ierr); 5719123193aSHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 5729123193aSHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 5739123193aSHong Zhang PetscFunctionReturn(0); 5749123193aSHong Zhang } 575cf1a0672SHong Zhang 576f8efd71cSHong Zhang PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C) 5771634b032SHong Zhang { 5781634b032SHong Zhang PetscErrorCode ierr; 579cf1a0672SHong Zhang Mat_MPIAIJ *a = (Mat_MPIAIJ*)A->data,*c=(Mat_MPIAIJ*)C->data; 580cf1a0672SHong Zhang Mat_SeqAIJ *ad = (Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data; 581cf1a0672SHong Zhang Mat_SeqAIJ *cd = (Mat_SeqAIJ*)(c->A)->data,*co=(Mat_SeqAIJ*)(c->B)->data; 582cf1a0672SHong Zhang PetscInt *adi = ad->i,*adj,*aoi=ao->i,*aoj; 583cf1a0672SHong Zhang PetscScalar *ada,*aoa,*cda=cd->a,*coa=co->a; 584cf1a0672SHong Zhang Mat_SeqAIJ *p_loc,*p_oth; 585cf1a0672SHong Zhang PetscInt *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pj; 58629825010SHong Zhang PetscScalar *pa_loc,*pa_oth,*pa,valtmp,*ca; 58729825010SHong Zhang PetscInt cm = C->rmap->n,anz,pnz; 5883cdca5ebSHong Zhang Mat_APMPI *ptap = c->ap; 58980da8df7SHong Zhang PetscScalar *apa_sparse; 590cf1a0672SHong Zhang PetscInt *api,*apj,*apJ,i,j,k,row; 59129825010SHong Zhang PetscInt cstart = C->cmap->rstart; 59229825010SHong Zhang PetscInt cdnz,conz,k0,k1,nextp; 593a7c7454dSHong Zhang MPI_Comm comm; 594a7c7454dSHong Zhang PetscMPIInt size; 5951634b032SHong Zhang 5961634b032SHong Zhang PetscFunctionBegin; 59780da8df7SHong Zhang ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr); 598a7c7454dSHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 599a7c7454dSHong Zhang 60080da8df7SHong Zhang if (!ptap) { 60180da8df7SHong Zhang SETERRQ(comm,PETSC_ERR_ARG_WRONGSTATE,"AP cannot be reused. Do not call MatFreeIntermediateDataStructures() or use '-mat_freeintermediatedatastructures'"); 60280da8df7SHong Zhang } 60380da8df7SHong Zhang apa_sparse = ptap->apa; 60480da8df7SHong Zhang 605cf1a0672SHong Zhang /* 1) get P_oth = ptap->P_oth and P_loc = ptap->P_loc */ 606cf1a0672SHong Zhang /*-----------------------------------------------------*/ 607cf1a0672SHong Zhang /* update numerical values of P_oth and P_loc */ 608b7f45c76SHong Zhang ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 609cf1a0672SHong Zhang ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 610cf1a0672SHong Zhang 611cf1a0672SHong Zhang /* 2) compute numeric C_loc = A_loc*P = Ad*P_loc + Ao*P_oth */ 612cf1a0672SHong Zhang /*----------------------------------------------------------*/ 613cf1a0672SHong Zhang /* get data from symbolic products */ 614cf1a0672SHong Zhang p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 615cf1a0672SHong Zhang pi_loc = p_loc->i; pj_loc = p_loc->j; pa_loc = p_loc->a; 616a7c7454dSHong Zhang if (size >1) { 617a7c7454dSHong Zhang p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 618cf1a0672SHong Zhang pi_oth = p_oth->i; pj_oth = p_oth->j; pa_oth = p_oth->a; 61920e3dc75SHong Zhang } else { 62020e3dc75SHong Zhang p_oth = NULL; pi_oth = NULL; pj_oth = NULL; pa_oth = NULL; 621a7c7454dSHong Zhang } 622cf1a0672SHong Zhang 62329825010SHong Zhang api = ptap->api; 62429825010SHong Zhang apj = ptap->apj; 625cf1a0672SHong Zhang for (i=0; i<cm; i++) { 62629825010SHong Zhang apJ = apj + api[i]; 62729825010SHong Zhang 628cf1a0672SHong Zhang /* diagonal portion of A */ 629cf1a0672SHong Zhang anz = adi[i+1] - adi[i]; 630cf1a0672SHong Zhang adj = ad->j + adi[i]; 631cf1a0672SHong Zhang ada = ad->a + adi[i]; 632cf1a0672SHong Zhang for (j=0; j<anz; j++) { 633cf1a0672SHong Zhang row = adj[j]; 634cf1a0672SHong Zhang pnz = pi_loc[row+1] - pi_loc[row]; 635cf1a0672SHong Zhang pj = pj_loc + pi_loc[row]; 636cf1a0672SHong Zhang pa = pa_loc + pi_loc[row]; 63729825010SHong Zhang /* perform sparse axpy */ 638cf1a0672SHong Zhang valtmp = ada[j]; 63929825010SHong Zhang nextp = 0; 64029825010SHong Zhang for (k=0; nextp<pnz; k++) { 64129825010SHong Zhang if (apJ[k] == pj[nextp]) { /* column of AP == column of P */ 64229825010SHong Zhang apa_sparse[k] += valtmp*pa[nextp++]; 64329825010SHong Zhang } 6441634b032SHong Zhang } 645cf1a0672SHong Zhang ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr); 646cf1a0672SHong Zhang } 6471634b032SHong Zhang 648cf1a0672SHong Zhang /* off-diagonal portion of A */ 649cf1a0672SHong Zhang anz = aoi[i+1] - aoi[i]; 650cf1a0672SHong Zhang aoj = ao->j + aoi[i]; 651cf1a0672SHong Zhang aoa = ao->a + aoi[i]; 652cf1a0672SHong Zhang for (j=0; j<anz; j++) { 653cf1a0672SHong Zhang row = aoj[j]; 654cf1a0672SHong Zhang pnz = pi_oth[row+1] - pi_oth[row]; 655cf1a0672SHong Zhang pj = pj_oth + pi_oth[row]; 656cf1a0672SHong Zhang pa = pa_oth + pi_oth[row]; 65729825010SHong Zhang /* perform sparse axpy */ 658cf1a0672SHong Zhang valtmp = aoa[j]; 65929825010SHong Zhang nextp = 0; 66029825010SHong Zhang for (k=0; nextp<pnz; k++) { 66129825010SHong Zhang if (apJ[k] == pj[nextp]) { /* column of AP == column of P */ 66229825010SHong Zhang apa_sparse[k] += valtmp*pa[nextp++]; 66329825010SHong Zhang } 664cf1a0672SHong Zhang } 665cf1a0672SHong Zhang ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr); 666cf1a0672SHong Zhang } 6671634b032SHong Zhang 668cf1a0672SHong Zhang /* set values in C */ 669cf1a0672SHong Zhang cdnz = cd->i[i+1] - cd->i[i]; 670cf1a0672SHong Zhang conz = co->i[i+1] - co->i[i]; 6711634b032SHong Zhang 672cf1a0672SHong Zhang /* 1st off-diagoanl part of C */ 673cf1a0672SHong Zhang ca = coa + co->i[i]; 674cf1a0672SHong Zhang k = 0; 675cf1a0672SHong Zhang for (k0=0; k0<conz; k0++) { 676cf1a0672SHong Zhang if (apJ[k] >= cstart) break; 67729825010SHong Zhang ca[k0] = apa_sparse[k]; 67829825010SHong Zhang apa_sparse[k] = 0.0; 679cf1a0672SHong Zhang k++; 680cf1a0672SHong Zhang } 6811634b032SHong Zhang 682cf1a0672SHong Zhang /* diagonal part of C */ 683cf1a0672SHong Zhang ca = cda + cd->i[i]; 684cf1a0672SHong Zhang for (k1=0; k1<cdnz; k1++) { 68529825010SHong Zhang ca[k1] = apa_sparse[k]; 68629825010SHong Zhang apa_sparse[k] = 0.0; 687cf1a0672SHong Zhang k++; 688cf1a0672SHong Zhang } 689cf1a0672SHong Zhang 690cf1a0672SHong Zhang /* 2nd off-diagoanl part of C */ 691cf1a0672SHong Zhang ca = coa + co->i[i]; 692cf1a0672SHong Zhang for (; k0<conz; k0++) { 69329825010SHong Zhang ca[k0] = apa_sparse[k]; 69429825010SHong Zhang apa_sparse[k] = 0.0; 695cf1a0672SHong Zhang k++; 696cf1a0672SHong Zhang } 697cf1a0672SHong Zhang } 698cf1a0672SHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 699cf1a0672SHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 70080da8df7SHong Zhang 701*7d0a89b7SHong Zhang if (ptap->freestruct) { 70280da8df7SHong Zhang ierr = MatFreeIntermediateDataStructures(C);CHKERRQ(ierr); 70380da8df7SHong Zhang } 704cf1a0672SHong Zhang PetscFunctionReturn(0); 705cf1a0672SHong Zhang } 706cf1a0672SHong Zhang 7070fc8cf34SHong Zhang /* same as MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(), except using LLCondensed to avoid O(BN) memory requirement */ 708b2405163SHong Zhang PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C) 709cf1a0672SHong Zhang { 710cf1a0672SHong Zhang PetscErrorCode ierr; 711ce94432eSBarry Smith MPI_Comm comm; 712a7c7454dSHong Zhang PetscMPIInt size; 713cf1a0672SHong Zhang Mat Cmpi; 7143cdca5ebSHong Zhang Mat_APMPI *ptap; 7150298fd71SBarry Smith PetscFreeSpaceList free_space = NULL,current_space=NULL; 716cf1a0672SHong Zhang Mat_MPIAIJ *a = (Mat_MPIAIJ*)A->data,*c; 717cf1a0672SHong Zhang Mat_SeqAIJ *ad = (Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_loc,*p_oth; 7187dbaa2c6Sandi selinger PetscInt *pi_loc,*pj_loc,*pi_oth,*pj_oth,*dnz,*onz; 719cf1a0672SHong Zhang PetscInt *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,rstart=A->rmap->rstart; 7208a210135Sandi selinger PetscInt i,pnz,row,*api,*apj,*Jptr,apnz,nspacedouble=0,j,nzi,*lnk,apnz_max=0; 721e5ab8c0aSKarl Rupp PetscInt am=A->rmap->n,pn=P->cmap->n,pm=P->rmap->n,lsize=pn+20; 722cf1a0672SHong Zhang PetscReal afill; 723cf1a0672SHong Zhang PetscScalar *apa; 724b92f168fSBarry Smith MatType mtype; 725cf1a0672SHong Zhang 726cf1a0672SHong Zhang PetscFunctionBegin; 727ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 728a7c7454dSHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 729a7c7454dSHong Zhang 7303cdca5ebSHong Zhang /* create struct Mat_APMPI and attached it to C later */ 731b00a9115SJed Brown ierr = PetscNew(&ptap);CHKERRQ(ierr); 732cf1a0672SHong Zhang 733cf1a0672SHong Zhang /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */ 734b7f45c76SHong Zhang ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 73519f0e57aSHong Zhang 736cf1a0672SHong Zhang /* get P_loc by taking all local rows of P */ 737cf1a0672SHong Zhang ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 738cf1a0672SHong Zhang 739cf1a0672SHong Zhang p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 740cf1a0672SHong Zhang pi_loc = p_loc->i; pj_loc = p_loc->j; 741a7c7454dSHong Zhang if (size > 1) { 742a7c7454dSHong Zhang p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 74320e3dc75SHong Zhang pi_oth = p_oth->i; pj_oth = p_oth->j; 744968382fdSHong Zhang } else { 74520e3dc75SHong Zhang p_oth = NULL; 74620e3dc75SHong Zhang pi_oth = NULL; pj_oth = NULL; 747a7c7454dSHong Zhang } 748cf1a0672SHong Zhang 749cf1a0672SHong Zhang /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */ 750cf1a0672SHong Zhang /*-------------------------------------------------------------------*/ 751854ce69bSBarry Smith ierr = PetscMalloc1(am+2,&api);CHKERRQ(ierr); 752cf1a0672SHong Zhang ptap->api = api; 753cf1a0672SHong Zhang api[0] = 0; 754cf1a0672SHong Zhang 7558a210135Sandi selinger ierr = PetscLLCondensedCreate_Scalable(lsize,&lnk);CHKERRQ(ierr); 756cf1a0672SHong Zhang 757cf1a0672SHong Zhang /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */ 758f91af8c7SBarry Smith ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(adi[am],PetscIntSumTruncate(aoi[am],pi_loc[pm]))),&free_space);CHKERRQ(ierr); 759cf1a0672SHong Zhang current_space = free_space; 760cf1a0672SHong Zhang ierr = MatPreallocateInitialize(comm,am,pn,dnz,onz);CHKERRQ(ierr); 761cf1a0672SHong Zhang for (i=0; i<am; i++) { 762cf1a0672SHong Zhang /* diagonal portion of A */ 763cf1a0672SHong Zhang nzi = adi[i+1] - adi[i]; 764cf1a0672SHong Zhang for (j=0; j<nzi; j++) { 765cf1a0672SHong Zhang row = *adj++; 766cf1a0672SHong Zhang pnz = pi_loc[row+1] - pi_loc[row]; 767cf1a0672SHong Zhang Jptr = pj_loc + pi_loc[row]; 7688a210135Sandi selinger /* Expand list if it is not long enough */ 7698328dd82Sandi selinger if (pnz+apnz_max > lsize) { 7708328dd82Sandi selinger lsize = pnz+apnz_max; 77126465080Sandi selinger ierr = PetscLLCondensedExpand_Scalable(lsize, &lnk);CHKERRQ(ierr); 7728a210135Sandi selinger } 773cf1a0672SHong Zhang /* add non-zero cols of P into the sorted linked list lnk */ 774f84c536bSHong Zhang ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr); 7758a210135Sandi selinger apnz = *lnk; /* The first element in the list is the number of items in the list */ 7768a210135Sandi selinger api[i+1] = api[i] + apnz; 7778a210135Sandi selinger if (apnz > apnz_max) apnz_max = apnz; 778cf1a0672SHong Zhang } 779cf1a0672SHong Zhang /* off-diagonal portion of A */ 780cf1a0672SHong Zhang nzi = aoi[i+1] - aoi[i]; 781cf1a0672SHong Zhang for (j=0; j<nzi; j++) { 782cf1a0672SHong Zhang row = *aoj++; 783cf1a0672SHong Zhang pnz = pi_oth[row+1] - pi_oth[row]; 784cf1a0672SHong Zhang Jptr = pj_oth + pi_oth[row]; 7858a210135Sandi selinger /* Expand list if it is not long enough */ 7868328dd82Sandi selinger if (pnz+apnz_max > lsize) { 7878328dd82Sandi selinger lsize = pnz + apnz_max; 78826465080Sandi selinger ierr = PetscLLCondensedExpand_Scalable(lsize, &lnk);CHKERRQ(ierr); 7898a210135Sandi selinger } 7908a210135Sandi selinger /* add non-zero cols of P into the sorted linked list lnk */ 791f84c536bSHong Zhang ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr); 7928a210135Sandi selinger apnz = *lnk; /* The first element in the list is the number of items in the list */ 7938a210135Sandi selinger api[i+1] = api[i] + apnz; 7948a210135Sandi selinger if (apnz > apnz_max) apnz_max = apnz; 795cf1a0672SHong Zhang } 796f84c536bSHong Zhang apnz = *lnk; 797cf1a0672SHong Zhang api[i+1] = api[i] + apnz; 7988a210135Sandi selinger if (apnz > apnz_max) apnz_max = apnz; 799cf1a0672SHong Zhang 800cf1a0672SHong Zhang /* if free space is not available, double the total space in the list */ 801cf1a0672SHong Zhang if (current_space->local_remaining<apnz) { 802f91af8c7SBarry Smith ierr = PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),¤t_space);CHKERRQ(ierr); 803cf1a0672SHong Zhang nspacedouble++; 804cf1a0672SHong Zhang } 805cf1a0672SHong Zhang 806cf1a0672SHong Zhang /* Copy data into free space, then initialize lnk */ 807f84c536bSHong Zhang ierr = PetscLLCondensedClean_Scalable(apnz,current_space->array,lnk);CHKERRQ(ierr); 808cf1a0672SHong Zhang ierr = MatPreallocateSet(i+rstart,apnz,current_space->array,dnz,onz);CHKERRQ(ierr); 8092205254eSKarl Rupp 810cf1a0672SHong Zhang current_space->array += apnz; 811cf1a0672SHong Zhang current_space->local_used += apnz; 812cf1a0672SHong Zhang current_space->local_remaining -= apnz; 813cf1a0672SHong Zhang } 814cf1a0672SHong Zhang 815cf1a0672SHong Zhang /* Allocate space for apj, initialize apj, and */ 816cf1a0672SHong Zhang /* destroy list of free space and other temporary array(s) */ 817854ce69bSBarry Smith ierr = PetscMalloc1(api[am]+1,&ptap->apj);CHKERRQ(ierr); 818cf1a0672SHong Zhang apj = ptap->apj; 819cf1a0672SHong Zhang ierr = PetscFreeSpaceContiguous(&free_space,ptap->apj);CHKERRQ(ierr); 820f84c536bSHong Zhang ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr); 821cf1a0672SHong Zhang 822cf1a0672SHong Zhang /* create and assemble symbolic parallel matrix Cmpi */ 823cf1a0672SHong Zhang /*----------------------------------------------------*/ 824cf1a0672SHong Zhang ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr); 825cf1a0672SHong Zhang ierr = MatSetSizes(Cmpi,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 82633d57670SJed Brown ierr = MatSetBlockSizesFromMats(Cmpi,A,P);CHKERRQ(ierr); 827b92f168fSBarry Smith ierr = MatGetType(A,&mtype);CHKERRQ(ierr); 828b92f168fSBarry Smith ierr = MatSetType(Cmpi,mtype);CHKERRQ(ierr); 829cf1a0672SHong Zhang ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr); 830cf1a0672SHong Zhang 83129825010SHong Zhang /* malloc apa for assembly Cmpi */ 8321795a4d1SJed Brown ierr = PetscCalloc1(apnz_max,&apa);CHKERRQ(ierr); 83329825010SHong Zhang ptap->apa = apa; 83417f75172Sandi selinger 835904d1e70Sandi selinger ierr = MatSetValues_MPIAIJ_CopyFromCSRFormat_Symbolic(Cmpi, apj, api);CHKERRQ(ierr); 836cf1a0672SHong Zhang ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 837cf1a0672SHong Zhang ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 8387dbaa2c6Sandi selinger ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 839cf1a0672SHong Zhang 840cf1a0672SHong Zhang ptap->destroy = Cmpi->ops->destroy; 841cf1a0672SHong Zhang ptap->duplicate = Cmpi->ops->duplicate; 842f8efd71cSHong Zhang Cmpi->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ; 843cf1a0672SHong Zhang Cmpi->ops->destroy = MatDestroy_MPIAIJ_MatMatMult; 8444624976aSHong Zhang Cmpi->ops->freeintermediatedatastructures = MatFreeIntermediateDataStructures_MPIAIJ_AP; 845cf1a0672SHong Zhang 846cf1a0672SHong Zhang /* attach the supporting struct to Cmpi for reuse */ 847cf1a0672SHong Zhang c = (Mat_MPIAIJ*)Cmpi->data; 8483cdca5ebSHong Zhang c->ap = ptap; 849cf1a0672SHong Zhang *C = Cmpi; 850cf1a0672SHong Zhang 851cf1a0672SHong Zhang /* set MatInfo */ 852118727c9SMark F. Adams afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1) + 1.e-5; 853cf1a0672SHong Zhang if (afill < 1.0) afill = 1.0; 854cf1a0672SHong Zhang Cmpi->info.mallocs = nspacedouble; 855cf1a0672SHong Zhang Cmpi->info.fill_ratio_given = fill; 856cf1a0672SHong Zhang Cmpi->info.fill_ratio_needed = afill; 857cf1a0672SHong Zhang 858cf1a0672SHong Zhang #if defined(PETSC_USE_INFO) 859cf1a0672SHong Zhang if (api[am]) { 86057622a8eSBarry Smith ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr); 86157622a8eSBarry Smith ierr = PetscInfo1(Cmpi,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr); 862cf1a0672SHong Zhang } else { 863cf1a0672SHong Zhang ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr); 864cf1a0672SHong Zhang } 865cf1a0672SHong Zhang #endif 8661634b032SHong Zhang PetscFunctionReturn(0); 8671634b032SHong Zhang } 868f96d37fcSHong Zhang 869a9d7e0ccSandi selinger /* This function is needed for the seqMPI matrix-matrix multiplication. */ 870a9d7e0ccSandi selinger /* Three input arrays are merged to one output array. The size of the */ 871a9d7e0ccSandi selinger /* output array is also output. Duplicate entries only show up once. */ 872a9d7e0ccSandi selinger static void Merge3SortedArrays(PetscInt size1, PetscInt *in1, 873a9d7e0ccSandi selinger PetscInt size2, PetscInt *in2, 874a9d7e0ccSandi selinger PetscInt size3, PetscInt *in3, 875a9d7e0ccSandi selinger PetscInt *size4, PetscInt *out) 876a9d7e0ccSandi selinger { 877a9d7e0ccSandi selinger int i = 0, j = 0, k = 0, l = 0; 878a9d7e0ccSandi selinger 879a9d7e0ccSandi selinger /* Traverse all three arrays */ 880a9d7e0ccSandi selinger while (i<size1 && j<size2 && k<size3) { 881a9d7e0ccSandi selinger if (in1[i] < in2[j] && in1[i] < in3[k]) { 882a9d7e0ccSandi selinger out[l++] = in1[i++]; 883a9d7e0ccSandi selinger } 884a9d7e0ccSandi selinger else if(in2[j] < in1[i] && in2[j] < in3[k]) { 885a9d7e0ccSandi selinger out[l++] = in2[j++]; 886a9d7e0ccSandi selinger } 887a9d7e0ccSandi selinger else if(in3[k] < in1[i] && in3[k] < in2[j]) { 888a9d7e0ccSandi selinger out[l++] = in3[k++]; 889a9d7e0ccSandi selinger } 890a9d7e0ccSandi selinger else if(in1[i] == in2[j] && in1[i] < in3[k]) { 891a9d7e0ccSandi selinger out[l++] = in1[i]; 892a9d7e0ccSandi selinger i++, j++; 893a9d7e0ccSandi selinger } 894a9d7e0ccSandi selinger else if(in1[i] == in3[k] && in1[i] < in2[j]) { 895a9d7e0ccSandi selinger out[l++] = in1[i]; 896a9d7e0ccSandi selinger i++, k++; 897a9d7e0ccSandi selinger } 898a9d7e0ccSandi selinger else if(in3[k] == in2[j] && in2[j] < in1[i]) { 899a9d7e0ccSandi selinger out[l++] = in2[j]; 900a9d7e0ccSandi selinger k++, j++; 901a9d7e0ccSandi selinger } 902a9d7e0ccSandi selinger else if(in1[i] == in2[j] && in1[i] == in3[k]) { 903a9d7e0ccSandi selinger out[l++] = in1[i]; 904a9d7e0ccSandi selinger i++, j++, k++; 905a9d7e0ccSandi selinger } 906a9d7e0ccSandi selinger } 907a9d7e0ccSandi selinger 908a9d7e0ccSandi selinger /* Traverse two remaining arrays */ 909a9d7e0ccSandi selinger while (i<size1 && j<size2) { 910a9d7e0ccSandi selinger if (in1[i] < in2[j]) { 911a9d7e0ccSandi selinger out[l++] = in1[i++]; 912a9d7e0ccSandi selinger } 913a9d7e0ccSandi selinger else if(in1[i] > in2[j]) { 914a9d7e0ccSandi selinger out[l++] = in2[j++]; 915a9d7e0ccSandi selinger } 916a9d7e0ccSandi selinger else { 917a9d7e0ccSandi selinger out[l++] = in1[i]; 918a9d7e0ccSandi selinger i++, j++; 919a9d7e0ccSandi selinger } 920a9d7e0ccSandi selinger } 921a9d7e0ccSandi selinger 922a9d7e0ccSandi selinger while (i<size1 && k<size3) { 923a9d7e0ccSandi selinger if (in1[i] < in3[k]) { 924a9d7e0ccSandi selinger out[l++] = in1[i++]; 925a9d7e0ccSandi selinger } 926a9d7e0ccSandi selinger else if(in1[i] > in3[k]) { 927a9d7e0ccSandi selinger out[l++] = in3[k++]; 928a9d7e0ccSandi selinger } 929a9d7e0ccSandi selinger else { 930a9d7e0ccSandi selinger out[l++] = in1[i]; 931a9d7e0ccSandi selinger i++, k++; 932a9d7e0ccSandi selinger } 933a9d7e0ccSandi selinger } 934a9d7e0ccSandi selinger 935a9d7e0ccSandi selinger while (k<size3 && j<size2) { 936a9d7e0ccSandi selinger if (in3[k] < in2[j]) { 937a9d7e0ccSandi selinger out[l++] = in3[k++]; 938a9d7e0ccSandi selinger } 939a9d7e0ccSandi selinger else if(in3[k] > in2[j]) { 940a9d7e0ccSandi selinger out[l++] = in2[j++]; 941a9d7e0ccSandi selinger } 942a9d7e0ccSandi selinger else { 943a9d7e0ccSandi selinger out[l++] = in3[k]; 944a9d7e0ccSandi selinger k++, j++; 945a9d7e0ccSandi selinger } 946a9d7e0ccSandi selinger } 947a9d7e0ccSandi selinger 948a9d7e0ccSandi selinger /* Traverse one remaining array */ 949a9d7e0ccSandi selinger while (i<size1) out[l++] = in1[i++]; 950a9d7e0ccSandi selinger while (j<size2) out[l++] = in2[j++]; 951a9d7e0ccSandi selinger while (k<size3) out[l++] = in3[k++]; 952a9d7e0ccSandi selinger 953a9d7e0ccSandi selinger *size4 = l; 954a9d7e0ccSandi selinger } 955a9d7e0ccSandi selinger 95673b67375Sandi selinger /* This matrix-matrix multiplication algorithm divides the multiplication into three multiplications and */ 95773b67375Sandi selinger /* adds up the products. Two of these three multiplications are performed with existing (sequential) */ 95873b67375Sandi selinger /* matrix-matrix multiplications. */ 959a9d7e0ccSandi selinger #undef __FUNCT__ 960a9d7e0ccSandi selinger #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ_seqMPI" 961a9d7e0ccSandi selinger PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_seqMPI(Mat A, Mat P, PetscReal fill, Mat *C) 962a9d7e0ccSandi selinger { 963a9d7e0ccSandi selinger PetscErrorCode ierr; 964a9d7e0ccSandi selinger MPI_Comm comm; 965a9d7e0ccSandi selinger PetscMPIInt size; 966a9d7e0ccSandi selinger Mat Cmpi; 9673cdca5ebSHong Zhang Mat_APMPI *ptap; 968a9d7e0ccSandi selinger PetscFreeSpaceList free_space_diag=NULL, current_space=NULL; 969a9d7e0ccSandi selinger Mat_MPIAIJ *a =(Mat_MPIAIJ*)A->data; 970a9d7e0ccSandi selinger Mat_SeqAIJ *ad =(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_loc; 971a9d7e0ccSandi selinger Mat_MPIAIJ *p =(Mat_MPIAIJ*)P->data; 972a9d7e0ccSandi selinger Mat_MPIAIJ *c; 973a9d7e0ccSandi selinger Mat_SeqAIJ *adpd_seq, *p_off, *aopoth_seq; 974a9d7e0ccSandi selinger PetscInt adponz, adpdnz; 9757dbaa2c6Sandi selinger PetscInt *pi_loc,*dnz,*onz; 976a9d7e0ccSandi selinger PetscInt *adi=ad->i,*adj=ad->j,*aoi=ao->i,rstart=A->rmap->rstart; 977a9d7e0ccSandi selinger PetscInt *lnk,i, i1=0,pnz,row,*adpoi,*adpoj, *api, *adpoJ, *aopJ, *apJ,*Jptr, aopnz, nspacedouble=0,j,nzi, 978a9d7e0ccSandi selinger *apj,apnz, *adpdi, *adpdj, *adpdJ, *poff_i, *poff_j, *j_temp, *aopothi, *aopothj; 979a016958eSandi selinger PetscInt am=A->rmap->n,pN=P->cmap->N,pn=P->cmap->n,pm=P->rmap->n, p_colstart, p_colend; 980a9d7e0ccSandi selinger PetscBT lnkbt; 981a9d7e0ccSandi selinger PetscScalar *apa; 982a9d7e0ccSandi selinger PetscReal afill; 983a9d7e0ccSandi selinger PetscMPIInt rank; 984a9d7e0ccSandi selinger Mat adpd, aopoth; 985b92f168fSBarry Smith MatType mtype; 986a9d7e0ccSandi selinger 987a9d7e0ccSandi selinger PetscFunctionBegin; 988a9d7e0ccSandi selinger ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 989a9d7e0ccSandi selinger ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 99073b67375Sandi selinger ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 991a016958eSandi selinger ierr = MatGetOwnershipRangeColumn(P, &p_colstart, &p_colend); CHKERRQ(ierr); 992a9d7e0ccSandi selinger 9933cdca5ebSHong Zhang /* create struct Mat_APMPI and attached it to C later */ 994a9d7e0ccSandi selinger ierr = PetscNew(&ptap);CHKERRQ(ierr); 995a9d7e0ccSandi selinger 996a9d7e0ccSandi selinger /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */ 997a9d7e0ccSandi selinger ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 998a9d7e0ccSandi selinger 999a9d7e0ccSandi selinger /* get P_loc by taking all local rows of P */ 1000a9d7e0ccSandi selinger ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 1001a9d7e0ccSandi selinger 1002a9d7e0ccSandi selinger 1003a9d7e0ccSandi selinger p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 1004a9d7e0ccSandi selinger pi_loc = p_loc->i; 1005a9d7e0ccSandi selinger 1006a9d7e0ccSandi selinger /* Allocate memory for the i arrays of the matrices A*P, A_diag*P_off and A_offd * P */ 1007a9d7e0ccSandi selinger ierr = PetscMalloc1(am+2,&api);CHKERRQ(ierr); 1008a9d7e0ccSandi selinger ierr = PetscMalloc1(am+2,&adpoi);CHKERRQ(ierr); 1009a9d7e0ccSandi selinger 1010a9d7e0ccSandi selinger adpoi[0] = 0; 1011a9d7e0ccSandi selinger ptap->api = api; 1012a9d7e0ccSandi selinger api[0] = 0; 1013a9d7e0ccSandi selinger 1014a9d7e0ccSandi selinger /* create and initialize a linked list, will be used for both A_diag * P_loc_off and A_offd * P_oth */ 1015a9d7e0ccSandi selinger ierr = PetscLLCondensedCreate(pN,pN,&lnk,&lnkbt);CHKERRQ(ierr); 1016a9d7e0ccSandi selinger ierr = MatPreallocateInitialize(comm,am,pn,dnz,onz);CHKERRQ(ierr); 1017a9d7e0ccSandi selinger 1018a9d7e0ccSandi selinger /* Symbolic calc of A_loc_diag * P_loc_diag */ 101973b67375Sandi selinger ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(a->A, p->A, fill, &adpd);CHKERRQ(ierr); 1020a9d7e0ccSandi selinger adpd_seq = (Mat_SeqAIJ*)((adpd)->data); 1021a9d7e0ccSandi selinger adpdi = adpd_seq->i; adpdj = adpd_seq->j; 1022a9d7e0ccSandi selinger p_off = (Mat_SeqAIJ*)((p->B)->data); 1023a9d7e0ccSandi selinger poff_i = p_off->i; poff_j = p_off->j; 1024a9d7e0ccSandi selinger 102573b67375Sandi selinger /* j_temp stores indices of a result row before they are added to the linked list */ 1026a9d7e0ccSandi selinger ierr = PetscMalloc1(pN+2,&j_temp);CHKERRQ(ierr); 1027a9d7e0ccSandi selinger 1028a9d7e0ccSandi selinger 1029a9d7e0ccSandi selinger /* Symbolic calc of the A_diag * p_loc_off */ 1030a9d7e0ccSandi selinger /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */ 1031a9d7e0ccSandi selinger ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(adi[am],PetscIntSumTruncate(aoi[am],pi_loc[pm]))),&free_space_diag);CHKERRQ(ierr); 1032a9d7e0ccSandi selinger current_space = free_space_diag; 1033a9d7e0ccSandi selinger 1034a9d7e0ccSandi selinger for (i=0; i<am; i++) { 1035a9d7e0ccSandi selinger /* A_diag * P_loc_off */ 1036a9d7e0ccSandi selinger nzi = adi[i+1] - adi[i]; 1037a9d7e0ccSandi selinger for (j=0; j<nzi; j++) { 1038a9d7e0ccSandi selinger row = *adj++; 1039a9d7e0ccSandi selinger pnz = poff_i[row+1] - poff_i[row]; 1040a9d7e0ccSandi selinger Jptr = poff_j + poff_i[row]; 1041a9d7e0ccSandi selinger for(i1 = 0; i1 < pnz; i1++) { 1042a9d7e0ccSandi selinger j_temp[i1] = p->garray[Jptr[i1]]; 1043a9d7e0ccSandi selinger } 1044a9d7e0ccSandi selinger /* add non-zero cols of P into the sorted linked list lnk */ 1045a9d7e0ccSandi selinger ierr = PetscLLCondensedAddSorted(pnz,j_temp,lnk,lnkbt);CHKERRQ(ierr); 1046a9d7e0ccSandi selinger } 1047a9d7e0ccSandi selinger 1048a9d7e0ccSandi selinger adponz = lnk[0]; 1049a9d7e0ccSandi selinger adpoi[i+1] = adpoi[i] + adponz; 1050a9d7e0ccSandi selinger 1051a9d7e0ccSandi selinger /* if free space is not available, double the total space in the list */ 1052a9d7e0ccSandi selinger if (current_space->local_remaining<adponz) { 1053a9d7e0ccSandi selinger ierr = PetscFreeSpaceGet(PetscIntSumTruncate(adponz,current_space->total_array_size),¤t_space);CHKERRQ(ierr); 1054a9d7e0ccSandi selinger nspacedouble++; 1055a9d7e0ccSandi selinger } 1056a9d7e0ccSandi selinger 1057a9d7e0ccSandi selinger /* Copy data into free space, then initialize lnk */ 1058a9d7e0ccSandi selinger ierr = PetscLLCondensedClean(pN,adponz,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 1059a9d7e0ccSandi selinger 1060a9d7e0ccSandi selinger current_space->array += adponz; 1061a9d7e0ccSandi selinger current_space->local_used += adponz; 1062a9d7e0ccSandi selinger current_space->local_remaining -= adponz; 1063a9d7e0ccSandi selinger } 1064a9d7e0ccSandi selinger 1065a9d7e0ccSandi selinger /* Symbolic calc of A_off * P_oth */ 106673b67375Sandi selinger ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(a->B, ptap->P_oth, fill, &aopoth);CHKERRQ(ierr); 1067a9d7e0ccSandi selinger aopoth_seq = (Mat_SeqAIJ*)((aopoth)->data); 1068a9d7e0ccSandi selinger aopothi = aopoth_seq->i; aopothj = aopoth_seq->j; 1069a9d7e0ccSandi selinger 1070a9d7e0ccSandi selinger /* Allocate space for apj, adpj, aopj, ... */ 1071a9d7e0ccSandi selinger /* destroy lists of free space and other temporary array(s) */ 1072a9d7e0ccSandi selinger 1073a9d7e0ccSandi selinger ierr = PetscMalloc1(aopothi[am] + adpoi[am] + adpdi[am]+2, &ptap->apj);CHKERRQ(ierr); 1074a9d7e0ccSandi selinger ierr = PetscMalloc1(adpoi[am]+2, &adpoj);CHKERRQ(ierr); 1075a9d7e0ccSandi selinger 1076a9d7e0ccSandi selinger /* Copy from linked list to j-array */ 1077a9d7e0ccSandi selinger ierr = PetscFreeSpaceContiguous(&free_space_diag,adpoj);CHKERRQ(ierr); 1078a9d7e0ccSandi selinger ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 1079a9d7e0ccSandi selinger 1080a9d7e0ccSandi selinger adpoJ = adpoj; 1081a9d7e0ccSandi selinger adpdJ = adpdj; 1082a9d7e0ccSandi selinger aopJ = aopothj; 1083a9d7e0ccSandi selinger apj = ptap->apj; 1084a9d7e0ccSandi selinger apJ = apj; /* still empty */ 1085a9d7e0ccSandi selinger 1086a9d7e0ccSandi selinger /* Merge j-arrays of A_off * P, A_diag * P_loc_off, and */ 1087a9d7e0ccSandi selinger /* A_diag * P_loc_diag to get A*P */ 1088a9d7e0ccSandi selinger for (i = 0; i < am; i++) { 1089a9d7e0ccSandi selinger aopnz = aopothi[i+1] - aopothi[i]; 1090a9d7e0ccSandi selinger adponz = adpoi[i+1] - adpoi[i]; 1091a9d7e0ccSandi selinger adpdnz = adpdi[i+1] - adpdi[i]; 1092a9d7e0ccSandi selinger 1093a9d7e0ccSandi selinger /* Correct indices from A_diag*P_diag */ 1094a9d7e0ccSandi selinger for(i1 = 0; i1 < adpdnz; i1++) { 1095a016958eSandi selinger adpdJ[i1] += p_colstart; 1096a9d7e0ccSandi selinger } 1097a9d7e0ccSandi selinger /* Merge j-arrays of A_diag * P_loc_off and A_diag * P_loc_diag and A_off * P_oth */ 1098a9d7e0ccSandi selinger Merge3SortedArrays(adponz, adpoJ, adpdnz, adpdJ, aopnz, aopJ, &apnz, apJ); 1099a9d7e0ccSandi selinger ierr = MatPreallocateSet(i+rstart, apnz, apJ, dnz, onz); CHKERRQ(ierr); 1100a9d7e0ccSandi selinger 1101a9d7e0ccSandi selinger aopJ += aopnz; 1102a9d7e0ccSandi selinger adpoJ += adponz; 1103a9d7e0ccSandi selinger adpdJ += adpdnz; 1104a9d7e0ccSandi selinger apJ += apnz; 1105a9d7e0ccSandi selinger api[i+1] = api[i] + apnz; 1106a9d7e0ccSandi selinger } 1107a9d7e0ccSandi selinger 1108a9d7e0ccSandi selinger /* malloc apa to store dense row A[i,:]*P */ 1109a9d7e0ccSandi selinger ierr = PetscCalloc1(pN+2,&apa);CHKERRQ(ierr); 1110a9d7e0ccSandi selinger 1111a9d7e0ccSandi selinger ptap->apa = apa; 1112a9d7e0ccSandi selinger /* create and assemble symbolic parallel matrix Cmpi */ 1113a9d7e0ccSandi selinger ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr); 1114a9d7e0ccSandi selinger ierr = MatSetSizes(Cmpi,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 1115a9d7e0ccSandi selinger ierr = MatSetBlockSizesFromMats(Cmpi,A,P);CHKERRQ(ierr); 1116b92f168fSBarry Smith ierr = MatGetType(A,&mtype);CHKERRQ(ierr); 1117b92f168fSBarry Smith ierr = MatSetType(Cmpi,mtype);CHKERRQ(ierr); 1118a9d7e0ccSandi selinger ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr); 111917f75172Sandi selinger 112017f75172Sandi selinger 1121904d1e70Sandi selinger ierr = MatSetValues_MPIAIJ_CopyFromCSRFormat_Symbolic(Cmpi, apj, api);CHKERRQ(ierr); 1122a9d7e0ccSandi selinger ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1123a9d7e0ccSandi selinger ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 11247dbaa2c6Sandi selinger ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 1125a9d7e0ccSandi selinger 112617f75172Sandi selinger 1127a9d7e0ccSandi selinger ptap->destroy = Cmpi->ops->destroy; 1128a9d7e0ccSandi selinger ptap->duplicate = Cmpi->ops->duplicate; 1129a9d7e0ccSandi selinger Cmpi->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable; 1130a9d7e0ccSandi selinger Cmpi->ops->destroy = MatDestroy_MPIAIJ_MatMatMult; 1131a9d7e0ccSandi selinger 1132a9d7e0ccSandi selinger /* attach the supporting struct to Cmpi for reuse */ 1133a9d7e0ccSandi selinger c = (Mat_MPIAIJ*)Cmpi->data; 11343cdca5ebSHong Zhang c->ap = ptap; 1135a9d7e0ccSandi selinger *C = Cmpi; 1136a9d7e0ccSandi selinger 1137a9d7e0ccSandi selinger /* set MatInfo */ 113873b67375Sandi selinger afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1) + 1.e-5; 113973b67375Sandi selinger if (afill < 1.0) afill = 1.0; 1140a9d7e0ccSandi selinger Cmpi->info.mallocs = nspacedouble; 1141a9d7e0ccSandi selinger Cmpi->info.fill_ratio_given = fill; 1142a9d7e0ccSandi selinger Cmpi->info.fill_ratio_needed = afill; 1143a9d7e0ccSandi selinger 1144a9d7e0ccSandi selinger #if defined(PETSC_USE_INFO) 1145a9d7e0ccSandi selinger if (api[am]) { 1146a9d7e0ccSandi selinger ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr); 1147a9d7e0ccSandi selinger ierr = PetscInfo1(Cmpi,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr); 1148a9d7e0ccSandi selinger } else { 1149a9d7e0ccSandi selinger ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr); 1150a9d7e0ccSandi selinger } 1151a9d7e0ccSandi selinger #endif 1152a9d7e0ccSandi selinger 115373b67375Sandi selinger ierr = MatDestroy(&aopoth);CHKERRQ(ierr); 115473b67375Sandi selinger ierr = MatDestroy(&adpd);CHKERRQ(ierr); 115573b67375Sandi selinger ierr = PetscFree(j_temp);CHKERRQ(ierr); 115673b67375Sandi selinger ierr = PetscFree(adpoj);CHKERRQ(ierr); 115773b67375Sandi selinger ierr = PetscFree(adpoi);CHKERRQ(ierr); 1158a9d7e0ccSandi selinger PetscFunctionReturn(0); 1159a9d7e0ccSandi selinger } 1160a9d7e0ccSandi selinger 1161a9d7e0ccSandi selinger 1162f96d37fcSHong Zhang /*-------------------------------------------------------------------------*/ 1163c5af039cSHong Zhang PetscErrorCode MatTransposeMatMult_MPIAIJ_MPIAIJ(Mat P,Mat A,MatReuse scall,PetscReal fill,Mat *C) 1164f96d37fcSHong Zhang { 1165f96d37fcSHong Zhang PetscErrorCode ierr; 1166c216dbf3SHong Zhang const char *algTypes[3] = {"scalable","nonscalable","matmatmult"}; 11674ede91e3SHong Zhang PetscInt aN=A->cmap->N,alg=1; /* set default algorithm */ 11684ede91e3SHong Zhang PetscBool flg; 1169f96d37fcSHong Zhang 1170f96d37fcSHong Zhang PetscFunctionBegin; 1171c216dbf3SHong Zhang if (scall == MAT_INITIAL_MATRIX) { 1172715a5346SHong Zhang ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"MatTransposeMatMult","Mat");CHKERRQ(ierr); 11734ede91e3SHong Zhang ierr = PetscOptionsEList("-mattransposematmult_via","Algorithmic approach","MatTransposeMatMult",algTypes,3,algTypes[1],&alg,&flg);CHKERRQ(ierr); 1174c216dbf3SHong Zhang ierr = PetscOptionsEnd();CHKERRQ(ierr); 1175c216dbf3SHong Zhang 1176669e9c86SHong Zhang ierr = PetscLogEventBegin(MAT_TransposeMatMultSymbolic,P,A,0,0);CHKERRQ(ierr); 1177c216dbf3SHong Zhang switch (alg) { 1178c216dbf3SHong Zhang case 1: 11794ede91e3SHong Zhang if (!flg && aN > 100000) { /* may switch to scalable algorithm as default */ 11804ede91e3SHong Zhang MatInfo Ainfo,Pinfo; 11814ede91e3SHong Zhang PetscInt nz_local; 11824ede91e3SHong Zhang PetscBool alg_scalable_loc=PETSC_FALSE,alg_scalable; 11834ede91e3SHong Zhang MPI_Comm comm; 11844ede91e3SHong Zhang 11854ede91e3SHong Zhang ierr = MatGetInfo(A,MAT_LOCAL,&Ainfo);CHKERRQ(ierr); 11864ede91e3SHong Zhang ierr = MatGetInfo(P,MAT_LOCAL,&Pinfo);CHKERRQ(ierr); 11874ede91e3SHong Zhang nz_local = (PetscInt)(Ainfo.nz_allocated + Pinfo.nz_allocated); /* estimated local nonzero entries */ 11884ede91e3SHong Zhang 11894ede91e3SHong Zhang if (aN > fill*nz_local) alg_scalable_loc = PETSC_TRUE; 11904ede91e3SHong Zhang ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 11914ede91e3SHong Zhang ierr = MPIU_Allreduce(&alg_scalable_loc,&alg_scalable,1,MPIU_BOOL,MPI_LOR,comm);CHKERRQ(ierr); 11924ede91e3SHong Zhang 11934ede91e3SHong Zhang if (alg_scalable) { 11944ede91e3SHong Zhang alg = 0; /* scalable algorithm would slower than nonscalable algorithm */ 11954ede91e3SHong Zhang ierr = MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(P,A,fill,C);CHKERRQ(ierr); 11964ede91e3SHong Zhang break; 11974ede91e3SHong Zhang } 11984ede91e3SHong Zhang } 11992bbb1c24SHong Zhang ierr = MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(P,A,fill,C);CHKERRQ(ierr); 1200c216dbf3SHong Zhang break; 1201c216dbf3SHong Zhang case 2: 1202275476c6SMatthew G. Knepley { 1203ab1b013aSHong Zhang Mat Pt; 12043cdca5ebSHong Zhang Mat_APMPI *ptap; 1205ab1b013aSHong Zhang Mat_MPIAIJ *c; 1206ab1b013aSHong Zhang ierr = MatTranspose(P,MAT_INITIAL_MATRIX,&Pt);CHKERRQ(ierr); 1207ab1b013aSHong Zhang ierr = MatMatMult(Pt,A,MAT_INITIAL_MATRIX,fill,C);CHKERRQ(ierr); 1208ab1b013aSHong Zhang c = (Mat_MPIAIJ*)(*C)->data; 12093cdca5ebSHong Zhang ptap = c->ap; 1210*7d0a89b7SHong Zhang if (ptap) { 1211ab1b013aSHong Zhang ptap->Pt = Pt; 121260552ceaSHong Zhang (*C)->ops->freeintermediatedatastructures = MatFreeIntermediateDataStructures_MPIAIJ_AP; 1213*7d0a89b7SHong Zhang } 1214*7d0a89b7SHong Zhang (*C)->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult; 1215c216dbf3SHong Zhang PetscFunctionReturn(0); 1216275476c6SMatthew G. Knepley } 1217c216dbf3SHong Zhang break; 12184ede91e3SHong Zhang default: /* scalable algorithm */ 12196da69ca6SHong Zhang ierr = MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(P,A,fill,C);CHKERRQ(ierr); 1220c216dbf3SHong Zhang break; 1221c216dbf3SHong Zhang } 1222669e9c86SHong Zhang ierr = PetscLogEventEnd(MAT_TransposeMatMultSymbolic,P,A,0,0);CHKERRQ(ierr); 1223*7d0a89b7SHong Zhang 1224*7d0a89b7SHong Zhang { 1225*7d0a89b7SHong Zhang Mat_MPIAIJ *c = (Mat_MPIAIJ*)(*C)->data; 1226*7d0a89b7SHong Zhang Mat_APMPI *ap = c->ap; 1227*7d0a89b7SHong Zhang ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)(*C)),((PetscObject)(*C))->prefix,"MatFreeIntermediateDataStructures","Mat");CHKERRQ(ierr); 1228*7d0a89b7SHong Zhang ap->freestruct = PETSC_FALSE; 1229*7d0a89b7SHong Zhang ierr = PetscOptionsBool("-mat_freeintermediatedatastructures","Free intermediate data structures", "MatFreeIntermediateDataStructures",ap->freestruct,&ap->freestruct, NULL);CHKERRQ(ierr); 1230*7d0a89b7SHong Zhang ierr = PetscOptionsEnd();CHKERRQ(ierr); 1231c216dbf3SHong Zhang } 1232*7d0a89b7SHong Zhang } 1233*7d0a89b7SHong Zhang 1234669e9c86SHong Zhang ierr = PetscLogEventBegin(MAT_TransposeMatMultNumeric,P,A,0,0);CHKERRQ(ierr); 1235c216dbf3SHong Zhang ierr = (*(*C)->ops->mattransposemultnumeric)(P,A,*C);CHKERRQ(ierr); 1236669e9c86SHong Zhang ierr = PetscLogEventEnd(MAT_TransposeMatMultNumeric,P,A,0,0);CHKERRQ(ierr); 1237ab1b013aSHong Zhang PetscFunctionReturn(0); 1238ab1b013aSHong Zhang } 1239ab1b013aSHong Zhang 1240c216dbf3SHong Zhang /* This routine only works when scall=MAT_REUSE_MATRIX! */ 1241c216dbf3SHong Zhang PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult(Mat P,Mat A,Mat C) 1242c216dbf3SHong Zhang { 1243c216dbf3SHong Zhang PetscErrorCode ierr; 12442bbb1c24SHong Zhang Mat_MPIAIJ *c=(Mat_MPIAIJ*)C->data; 12453cdca5ebSHong Zhang Mat_APMPI *ptap= c->ap; 124660552ceaSHong Zhang Mat Pt; 1247c216dbf3SHong Zhang 1248c216dbf3SHong Zhang PetscFunctionBegin; 124960552ceaSHong Zhang if (!ptap) { 125060552ceaSHong Zhang MPI_Comm comm; 125160552ceaSHong Zhang ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr); 125260552ceaSHong Zhang SETERRQ(comm,PETSC_ERR_ARG_WRONGSTATE,"PtA cannot be reused. Do not call MatFreeIntermediateDataStructures() or use '-mat_freeintermediatedatastructures'"); 125360552ceaSHong Zhang } 125460552ceaSHong Zhang 125560552ceaSHong Zhang Pt=ptap->Pt; 1256c216dbf3SHong Zhang ierr = MatTranspose(P,MAT_REUSE_MATRIX,&Pt);CHKERRQ(ierr); 1257c216dbf3SHong Zhang ierr = MatMatMultNumeric(Pt,A,C);CHKERRQ(ierr); 125860552ceaSHong Zhang 125960552ceaSHong Zhang /* supporting struct ptap consumes almost same amount of memory as C=PtAP, release it if C will not be updated by A and P */ 1260*7d0a89b7SHong Zhang if (ptap->freestruct) { 126160552ceaSHong Zhang ierr = MatFreeIntermediateDataStructures(C);CHKERRQ(ierr); 126260552ceaSHong Zhang } 1263f96d37fcSHong Zhang PetscFunctionReturn(0); 1264f96d37fcSHong Zhang } 1265f96d37fcSHong Zhang 1266c5af039cSHong Zhang /* This routine is modified from MatPtAPSymbolic_MPIAIJ_MPIAIJ() */ 1267912837bbSHong Zhang PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat P,Mat A,PetscReal fill,Mat *C) 1268912837bbSHong Zhang { 1269912837bbSHong Zhang PetscErrorCode ierr; 12703cdca5ebSHong Zhang Mat_APMPI *ptap; 12717f3e1f30SHong Zhang Mat_MPIAIJ *p=(Mat_MPIAIJ*)P->data,*c; 1272912837bbSHong Zhang MPI_Comm comm; 1273912837bbSHong Zhang PetscMPIInt size,rank; 1274912837bbSHong Zhang Mat Cmpi; 1275912837bbSHong Zhang PetscFreeSpaceList free_space=NULL,current_space=NULL; 12764ede91e3SHong Zhang PetscInt pn=P->cmap->n,aN=A->cmap->N,an=A->cmap->n; 12777f3e1f30SHong Zhang PetscInt *lnk,i,k,nsend; 1278912837bbSHong Zhang PetscBT lnkbt; 1279912837bbSHong Zhang PetscMPIInt tagi,tagj,*len_si,*len_s,*len_ri,icompleted=0,nrecv; 1280912837bbSHong Zhang PetscInt **buf_rj,**buf_ri,**buf_ri_k; 12817f3e1f30SHong Zhang PetscInt len,proc,*dnz,*onz,*owners,nzi; 1282912837bbSHong Zhang PetscInt nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci; 1283912837bbSHong Zhang MPI_Request *swaits,*rwaits; 1284912837bbSHong Zhang MPI_Status *sstatus,rstatus; 1285912837bbSHong Zhang PetscLayout rowmap; 1286912837bbSHong Zhang PetscInt *owners_co,*coi,*coj; /* i and j array of (p->B)^T*A*P - used in the communication */ 1287912837bbSHong Zhang PetscMPIInt *len_r,*id_r; /* array of length of comm->size, store send/recv matrix values */ 12887f3e1f30SHong Zhang PetscInt *Jptr,*prmap=p->garray,con,j,Crmax; 12894ede91e3SHong Zhang Mat_SeqAIJ *a_loc,*c_loc,*c_oth; 1290912837bbSHong Zhang PetscTable ta; 1291b92f168fSBarry Smith MatType mtype; 1292912837bbSHong Zhang 1293912837bbSHong Zhang PetscFunctionBegin; 1294912837bbSHong Zhang ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 1295912837bbSHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1296912837bbSHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 1297912837bbSHong Zhang 1298912837bbSHong Zhang /* create symbolic parallel matrix Cmpi */ 1299912837bbSHong Zhang ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr); 1300b92f168fSBarry Smith ierr = MatGetType(A,&mtype);CHKERRQ(ierr); 1301b92f168fSBarry Smith ierr = MatSetType(Cmpi,mtype);CHKERRQ(ierr); 1302912837bbSHong Zhang 1303912837bbSHong Zhang Cmpi->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable; 1304912837bbSHong Zhang 13053cdca5ebSHong Zhang /* create struct Mat_APMPI and attached it to C later */ 1306912837bbSHong Zhang ierr = PetscNew(&ptap);CHKERRQ(ierr); 1307912837bbSHong Zhang ptap->reuse = MAT_INITIAL_MATRIX; 1308912837bbSHong Zhang 1309912837bbSHong Zhang /* (0) compute Rd = Pd^T, Ro = Po^T */ 1310912837bbSHong Zhang /* --------------------------------- */ 1311912837bbSHong Zhang ierr = MatTranspose_SeqAIJ(p->A,MAT_INITIAL_MATRIX,&ptap->Rd);CHKERRQ(ierr); 1312912837bbSHong Zhang ierr = MatTranspose_SeqAIJ(p->B,MAT_INITIAL_MATRIX,&ptap->Ro);CHKERRQ(ierr); 1313912837bbSHong Zhang 1314669e9c86SHong Zhang /* (1) compute symbolic A_loc */ 1315669e9c86SHong Zhang /* ---------------------------*/ 13164ede91e3SHong Zhang ierr = MatMPIAIJGetLocalMat(A,MAT_INITIAL_MATRIX,&ptap->A_loc);CHKERRQ(ierr); 1317912837bbSHong Zhang 13184ede91e3SHong Zhang /* (2-1) compute symbolic C_oth = Ro*A_loc */ 1319912837bbSHong Zhang /* ------------------------------------ */ 13204ede91e3SHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Ro,ptap->A_loc,fill,&ptap->C_oth);CHKERRQ(ierr); 1321912837bbSHong Zhang 1322912837bbSHong Zhang /* (3) send coj of C_oth to other processors */ 1323912837bbSHong Zhang /* ------------------------------------------ */ 1324912837bbSHong Zhang /* determine row ownership */ 1325912837bbSHong Zhang ierr = PetscLayoutCreate(comm,&rowmap);CHKERRQ(ierr); 1326912837bbSHong Zhang rowmap->n = pn; 1327912837bbSHong Zhang rowmap->bs = 1; 1328912837bbSHong Zhang ierr = PetscLayoutSetUp(rowmap);CHKERRQ(ierr); 1329912837bbSHong Zhang owners = rowmap->range; 1330912837bbSHong Zhang 1331912837bbSHong Zhang /* determine the number of messages to send, their lengths */ 1332912837bbSHong Zhang ierr = PetscMalloc4(size,&len_s,size,&len_si,size,&sstatus,size+2,&owners_co);CHKERRQ(ierr); 1333912837bbSHong Zhang ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 1334912837bbSHong Zhang ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 1335912837bbSHong Zhang 1336912837bbSHong Zhang c_oth = (Mat_SeqAIJ*)ptap->C_oth->data; 1337912837bbSHong Zhang coi = c_oth->i; coj = c_oth->j; 1338912837bbSHong Zhang con = ptap->C_oth->rmap->n; 1339912837bbSHong Zhang proc = 0; 1340912837bbSHong Zhang for (i=0; i<con; i++) { 1341912837bbSHong Zhang while (prmap[i] >= owners[proc+1]) proc++; 13427f3e1f30SHong Zhang len_si[proc]++; /* num of rows in Co(=Pt*A) to be sent to [proc] */ 1343912837bbSHong Zhang len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */ 1344912837bbSHong Zhang } 1345912837bbSHong Zhang 1346912837bbSHong Zhang len = 0; /* max length of buf_si[], see (4) */ 1347912837bbSHong Zhang owners_co[0] = 0; 1348912837bbSHong Zhang nsend = 0; 1349912837bbSHong Zhang for (proc=0; proc<size; proc++) { 1350912837bbSHong Zhang owners_co[proc+1] = owners_co[proc] + len_si[proc]; 1351912837bbSHong Zhang if (len_s[proc]) { 1352912837bbSHong Zhang nsend++; 1353912837bbSHong Zhang len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */ 1354912837bbSHong Zhang len += len_si[proc]; 1355912837bbSHong Zhang } 1356912837bbSHong Zhang } 1357912837bbSHong Zhang 1358912837bbSHong Zhang /* determine the number and length of messages to receive for coi and coj */ 1359912837bbSHong Zhang ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&nrecv);CHKERRQ(ierr); 1360912837bbSHong Zhang ierr = PetscGatherMessageLengths2(comm,nsend,nrecv,len_s,len_si,&id_r,&len_r,&len_ri);CHKERRQ(ierr); 1361912837bbSHong Zhang 1362912837bbSHong Zhang /* post the Irecv and Isend of coj */ 1363912837bbSHong Zhang ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr); 1364912837bbSHong Zhang ierr = PetscPostIrecvInt(comm,tagj,nrecv,id_r,len_r,&buf_rj,&rwaits);CHKERRQ(ierr); 1365912837bbSHong Zhang ierr = PetscMalloc1(nsend+1,&swaits);CHKERRQ(ierr); 1366912837bbSHong Zhang for (proc=0, k=0; proc<size; proc++) { 1367912837bbSHong Zhang if (!len_s[proc]) continue; 1368912837bbSHong Zhang i = owners_co[proc]; 1369912837bbSHong Zhang ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr); 1370912837bbSHong Zhang k++; 1371912837bbSHong Zhang } 1372912837bbSHong Zhang 13737f3e1f30SHong Zhang /* (2-2) compute symbolic C_loc = Rd*A_loc */ 1374912837bbSHong Zhang /* ---------------------------------------- */ 13754ede91e3SHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Rd,ptap->A_loc,fill,&ptap->C_loc);CHKERRQ(ierr); 1376912837bbSHong Zhang c_loc = (Mat_SeqAIJ*)ptap->C_loc->data; 1377912837bbSHong Zhang 1378912837bbSHong Zhang /* receives coj are complete */ 1379912837bbSHong Zhang for (i=0; i<nrecv; i++) { 1380912837bbSHong Zhang ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 1381912837bbSHong Zhang } 1382912837bbSHong Zhang ierr = PetscFree(rwaits);CHKERRQ(ierr); 1383912837bbSHong Zhang if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);} 1384912837bbSHong Zhang 1385912837bbSHong Zhang /* add received column indices into ta to update Crmax */ 13864ede91e3SHong Zhang a_loc = (Mat_SeqAIJ*)(ptap->A_loc)->data; 1387669e9c86SHong Zhang 1388669e9c86SHong Zhang /* create and initialize a linked list */ 13894ede91e3SHong Zhang ierr = PetscTableCreate(an,aN,&ta);CHKERRQ(ierr); /* for compute Crmax */ 13904ede91e3SHong Zhang MatRowMergeMax_SeqAIJ(a_loc,ptap->A_loc->rmap->N,ta); 1391669e9c86SHong Zhang 1392912837bbSHong Zhang for (k=0; k<nrecv; k++) {/* k-th received message */ 1393912837bbSHong Zhang Jptr = buf_rj[k]; 1394912837bbSHong Zhang for (j=0; j<len_r[k]; j++) { 1395912837bbSHong Zhang ierr = PetscTableAdd(ta,*(Jptr+j)+1,1,INSERT_VALUES);CHKERRQ(ierr); 1396912837bbSHong Zhang } 1397912837bbSHong Zhang } 1398912837bbSHong Zhang ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr); 1399912837bbSHong Zhang ierr = PetscTableDestroy(&ta);CHKERRQ(ierr); 1400912837bbSHong Zhang 1401912837bbSHong Zhang /* (4) send and recv coi */ 1402912837bbSHong Zhang /*-----------------------*/ 1403912837bbSHong Zhang ierr = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr); 1404912837bbSHong Zhang ierr = PetscPostIrecvInt(comm,tagi,nrecv,id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr); 1405912837bbSHong Zhang ierr = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr); 1406912837bbSHong Zhang buf_si = buf_s; /* points to the beginning of k-th msg to be sent */ 1407912837bbSHong Zhang for (proc=0,k=0; proc<size; proc++) { 1408912837bbSHong Zhang if (!len_s[proc]) continue; 1409912837bbSHong Zhang /* form outgoing message for i-structure: 1410912837bbSHong Zhang buf_si[0]: nrows to be sent 1411912837bbSHong Zhang [1:nrows]: row index (global) 1412912837bbSHong Zhang [nrows+1:2*nrows+1]: i-structure index 1413912837bbSHong Zhang */ 1414912837bbSHong Zhang /*-------------------------------------------*/ 1415912837bbSHong Zhang nrows = len_si[proc]/2 - 1; /* num of rows in Co to be sent to [proc] */ 1416912837bbSHong Zhang buf_si_i = buf_si + nrows+1; 1417912837bbSHong Zhang buf_si[0] = nrows; 1418912837bbSHong Zhang buf_si_i[0] = 0; 1419912837bbSHong Zhang nrows = 0; 1420912837bbSHong Zhang for (i=owners_co[proc]; i<owners_co[proc+1]; i++) { 1421912837bbSHong Zhang nzi = coi[i+1] - coi[i]; 1422912837bbSHong Zhang buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */ 1423912837bbSHong Zhang buf_si[nrows+1] = prmap[i] -owners[proc]; /* local row index */ 1424912837bbSHong Zhang nrows++; 1425912837bbSHong Zhang } 1426912837bbSHong Zhang ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr); 1427912837bbSHong Zhang k++; 1428912837bbSHong Zhang buf_si += len_si[proc]; 1429912837bbSHong Zhang } 1430912837bbSHong Zhang for (i=0; i<nrecv; i++) { 1431912837bbSHong Zhang ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 1432912837bbSHong Zhang } 1433912837bbSHong Zhang ierr = PetscFree(rwaits);CHKERRQ(ierr); 1434912837bbSHong Zhang if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);} 1435912837bbSHong Zhang 1436912837bbSHong Zhang ierr = PetscFree4(len_s,len_si,sstatus,owners_co);CHKERRQ(ierr); 1437912837bbSHong Zhang ierr = PetscFree(len_ri);CHKERRQ(ierr); 1438912837bbSHong Zhang ierr = PetscFree(swaits);CHKERRQ(ierr); 1439912837bbSHong Zhang ierr = PetscFree(buf_s);CHKERRQ(ierr); 1440912837bbSHong Zhang 1441912837bbSHong Zhang /* (5) compute the local portion of Cmpi */ 1442912837bbSHong Zhang /* ------------------------------------------ */ 1443912837bbSHong Zhang /* set initial free space to be Crmax, sufficient for holding nozeros in each row of Cmpi */ 1444912837bbSHong Zhang ierr = PetscFreeSpaceGet(Crmax,&free_space);CHKERRQ(ierr); 1445912837bbSHong Zhang current_space = free_space; 1446912837bbSHong Zhang 1447912837bbSHong Zhang ierr = PetscMalloc3(nrecv,&buf_ri_k,nrecv,&nextrow,nrecv,&nextci);CHKERRQ(ierr); 1448912837bbSHong Zhang for (k=0; k<nrecv; k++) { 1449912837bbSHong Zhang buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 1450912837bbSHong Zhang nrows = *buf_ri_k[k]; 1451912837bbSHong Zhang nextrow[k] = buf_ri_k[k] + 1; /* next row number of k-th recved i-structure */ 1452912837bbSHong Zhang nextci[k] = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure */ 1453912837bbSHong Zhang } 1454912837bbSHong Zhang 14554ede91e3SHong Zhang ierr = MatPreallocateInitialize(comm,pn,an,dnz,onz);CHKERRQ(ierr); 14564ede91e3SHong Zhang ierr = PetscLLCondensedCreate(Crmax,aN,&lnk,&lnkbt);CHKERRQ(ierr); 1457912837bbSHong Zhang for (i=0; i<pn; i++) { 1458912837bbSHong Zhang /* add C_loc into Cmpi */ 1459912837bbSHong Zhang nzi = c_loc->i[i+1] - c_loc->i[i]; 1460912837bbSHong Zhang Jptr = c_loc->j + c_loc->i[i]; 1461912837bbSHong Zhang ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr); 1462912837bbSHong Zhang 1463912837bbSHong Zhang /* add received col data into lnk */ 1464912837bbSHong Zhang for (k=0; k<nrecv; k++) { /* k-th received message */ 1465912837bbSHong Zhang if (i == *nextrow[k]) { /* i-th row */ 1466912837bbSHong Zhang nzi = *(nextci[k]+1) - *nextci[k]; 1467912837bbSHong Zhang Jptr = buf_rj[k] + *nextci[k]; 1468912837bbSHong Zhang ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr); 1469912837bbSHong Zhang nextrow[k]++; nextci[k]++; 1470912837bbSHong Zhang } 1471912837bbSHong Zhang } 1472912837bbSHong Zhang nzi = lnk[0]; 1473912837bbSHong Zhang 1474912837bbSHong Zhang /* copy data into free space, then initialize lnk */ 14754ede91e3SHong Zhang ierr = PetscLLCondensedClean(aN,nzi,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 1476912837bbSHong Zhang ierr = MatPreallocateSet(i+owners[rank],nzi,current_space->array,dnz,onz);CHKERRQ(ierr); 1477912837bbSHong Zhang } 1478912837bbSHong Zhang ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr); 1479912837bbSHong Zhang ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 1480912837bbSHong Zhang ierr = PetscFreeSpaceDestroy(free_space);CHKERRQ(ierr); 1481912837bbSHong Zhang 1482912837bbSHong Zhang /* local sizes and preallocation */ 14834ede91e3SHong Zhang ierr = MatSetSizes(Cmpi,pn,an,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 1484912837bbSHong Zhang ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr); 1485912837bbSHong Zhang ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr); 1486912837bbSHong Zhang ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 1487912837bbSHong Zhang 1488912837bbSHong Zhang /* members in merge */ 1489912837bbSHong Zhang ierr = PetscFree(id_r);CHKERRQ(ierr); 1490912837bbSHong Zhang ierr = PetscFree(len_r);CHKERRQ(ierr); 1491912837bbSHong Zhang ierr = PetscFree(buf_ri[0]);CHKERRQ(ierr); 1492912837bbSHong Zhang ierr = PetscFree(buf_ri);CHKERRQ(ierr); 1493912837bbSHong Zhang ierr = PetscFree(buf_rj[0]);CHKERRQ(ierr); 1494912837bbSHong Zhang ierr = PetscFree(buf_rj);CHKERRQ(ierr); 1495912837bbSHong Zhang ierr = PetscLayoutDestroy(&rowmap);CHKERRQ(ierr); 1496912837bbSHong Zhang 1497912837bbSHong Zhang /* attach the supporting struct to Cmpi for reuse */ 1498912837bbSHong Zhang c = (Mat_MPIAIJ*)Cmpi->data; 14993cdca5ebSHong Zhang c->ap = ptap; 1500912837bbSHong Zhang ptap->destroy = Cmpi->ops->destroy; 1501912837bbSHong Zhang 1502912837bbSHong Zhang /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */ 1503912837bbSHong Zhang Cmpi->assembled = PETSC_FALSE; 1504912837bbSHong Zhang Cmpi->ops->destroy = MatDestroy_MPIAIJ_PtAP; 15054624976aSHong Zhang Cmpi->ops->freeintermediatedatastructures = MatFreeIntermediateDataStructures_MPIAIJ_AP; 15068a9b8717SHong Zhang 1507912837bbSHong Zhang *C = Cmpi; 1508912837bbSHong Zhang PetscFunctionReturn(0); 1509912837bbSHong Zhang } 1510912837bbSHong Zhang 1511912837bbSHong Zhang PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable(Mat P,Mat A,Mat C) 1512912837bbSHong Zhang { 1513912837bbSHong Zhang PetscErrorCode ierr; 1514669e9c86SHong Zhang Mat_MPIAIJ *p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data; 1515669e9c86SHong Zhang Mat_SeqAIJ *c_seq; 15163cdca5ebSHong Zhang Mat_APMPI *ptap = c->ap; 1517669e9c86SHong Zhang Mat A_loc,C_loc,C_oth; 1518912837bbSHong Zhang PetscInt i,rstart,rend,cm,ncols,row; 1519912837bbSHong Zhang const PetscInt *cols; 1520912837bbSHong Zhang const PetscScalar *vals; 1521912837bbSHong Zhang 1522912837bbSHong Zhang PetscFunctionBegin; 15238a9b8717SHong Zhang if (!ptap) { 15248a9b8717SHong Zhang MPI_Comm comm; 15258a9b8717SHong Zhang ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr); 15268a9b8717SHong Zhang SETERRQ(comm,PETSC_ERR_ARG_WRONGSTATE,"PtA cannot be reused. Do not call MatFreeIntermediateDataStructures() or use '-mat_freeintermediatedatastructures'"); 15278a9b8717SHong Zhang } 15288a9b8717SHong Zhang 1529912837bbSHong Zhang ierr = MatZeroEntries(C);CHKERRQ(ierr); 1530912837bbSHong Zhang 1531912837bbSHong Zhang if (ptap->reuse == MAT_REUSE_MATRIX) { 1532669e9c86SHong Zhang /* These matrices are obtained in MatTransposeMatMultSymbolic() */ 1533669e9c86SHong Zhang /* 1) get R = Pd^T, Ro = Po^T */ 1534669e9c86SHong Zhang /*----------------------------*/ 1535912837bbSHong Zhang ierr = MatTranspose_SeqAIJ(p->A,MAT_REUSE_MATRIX,&ptap->Rd);CHKERRQ(ierr); 1536912837bbSHong Zhang ierr = MatTranspose_SeqAIJ(p->B,MAT_REUSE_MATRIX,&ptap->Ro);CHKERRQ(ierr); 1537669e9c86SHong Zhang 1538669e9c86SHong Zhang /* 2) compute numeric A_loc */ 1539669e9c86SHong Zhang /*--------------------------*/ 15404ede91e3SHong Zhang ierr = MatMPIAIJGetLocalMat(A,MAT_REUSE_MATRIX,&ptap->A_loc);CHKERRQ(ierr); 1541912837bbSHong Zhang } 1542912837bbSHong Zhang 1543669e9c86SHong Zhang /* 3) C_loc = Rd*A_loc, C_oth = Ro*A_loc */ 15444ede91e3SHong Zhang A_loc = ptap->A_loc; 1545669e9c86SHong Zhang ierr = ((ptap->C_loc)->ops->matmultnumeric)(ptap->Rd,A_loc,ptap->C_loc);CHKERRQ(ierr); 1546669e9c86SHong Zhang ierr = ((ptap->C_oth)->ops->matmultnumeric)(ptap->Ro,A_loc,ptap->C_oth);CHKERRQ(ierr); 1547912837bbSHong Zhang C_loc = ptap->C_loc; 1548912837bbSHong Zhang C_oth = ptap->C_oth; 1549912837bbSHong Zhang 1550912837bbSHong Zhang /* add C_loc and Co to to C */ 1551912837bbSHong Zhang ierr = MatGetOwnershipRange(C,&rstart,&rend);CHKERRQ(ierr); 1552912837bbSHong Zhang 1553912837bbSHong Zhang /* C_loc -> C */ 1554912837bbSHong Zhang cm = C_loc->rmap->N; 1555912837bbSHong Zhang c_seq = (Mat_SeqAIJ*)C_loc->data; 1556912837bbSHong Zhang cols = c_seq->j; 1557912837bbSHong Zhang vals = c_seq->a; 1558912837bbSHong Zhang for (i=0; i<cm; i++) { 1559912837bbSHong Zhang ncols = c_seq->i[i+1] - c_seq->i[i]; 1560912837bbSHong Zhang row = rstart + i; 1561912837bbSHong Zhang ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr); 1562912837bbSHong Zhang cols += ncols; vals += ncols; 1563912837bbSHong Zhang } 1564912837bbSHong Zhang 1565912837bbSHong Zhang /* Co -> C, off-processor part */ 1566912837bbSHong Zhang cm = C_oth->rmap->N; 1567912837bbSHong Zhang c_seq = (Mat_SeqAIJ*)C_oth->data; 1568912837bbSHong Zhang cols = c_seq->j; 1569912837bbSHong Zhang vals = c_seq->a; 1570912837bbSHong Zhang for (i=0; i<cm; i++) { 1571912837bbSHong Zhang ncols = c_seq->i[i+1] - c_seq->i[i]; 1572912837bbSHong Zhang row = p->garray[i]; 1573912837bbSHong Zhang ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr); 1574912837bbSHong Zhang cols += ncols; vals += ncols; 1575912837bbSHong Zhang } 1576912837bbSHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1577912837bbSHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1578912837bbSHong Zhang 1579912837bbSHong Zhang ptap->reuse = MAT_REUSE_MATRIX; 15808a9b8717SHong Zhang 15818a9b8717SHong Zhang /* supporting struct ptap consumes almost same amount of memory as C=PtAP, release it if C will not be updated by A and P */ 1582*7d0a89b7SHong Zhang if (ptap->freestruct) { 15838a9b8717SHong Zhang ierr = MatFreeIntermediateDataStructures(C);CHKERRQ(ierr); 15848a9b8717SHong Zhang } 1585912837bbSHong Zhang PetscFunctionReturn(0); 1586912837bbSHong Zhang } 1587912837bbSHong Zhang 15886da69ca6SHong Zhang PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ(Mat P,Mat A,Mat C) 1589d6ab1dc0SHong Zhang { 1590d6ab1dc0SHong Zhang PetscErrorCode ierr; 1591d6ab1dc0SHong Zhang Mat_Merge_SeqsToMPI *merge; 1592d6ab1dc0SHong Zhang Mat_MPIAIJ *p =(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data; 1593d6ab1dc0SHong Zhang Mat_SeqAIJ *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data; 15943cdca5ebSHong Zhang Mat_APMPI *ptap; 1595e745005bSHong Zhang PetscInt *adj; 1596e745005bSHong Zhang PetscInt i,j,k,anz,pnz,row,*cj,nexta; 1597e745005bSHong Zhang MatScalar *ada,*ca,valtmp; 1598d6ab1dc0SHong Zhang PetscInt am =A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n; 1599ce94432eSBarry Smith MPI_Comm comm; 1600d6ab1dc0SHong Zhang PetscMPIInt size,rank,taga,*len_s; 1601d6ab1dc0SHong Zhang PetscInt *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci; 1602d6ab1dc0SHong Zhang PetscInt **buf_ri,**buf_rj; 1603d6ab1dc0SHong Zhang PetscInt cnz=0,*bj_i,*bi,*bj,bnz,nextcj; /* bi,bj,ba: local array of C(mpi mat) */ 1604d6ab1dc0SHong Zhang MPI_Request *s_waits,*r_waits; 1605d6ab1dc0SHong Zhang MPI_Status *status; 1606d6ab1dc0SHong Zhang MatScalar **abuf_r,*ba_i,*pA,*coa,*ba; 16074e201bd7SHong Zhang PetscInt *ai,*aj,*coi,*coj,*poJ,*pdJ; 1608d6ab1dc0SHong Zhang Mat A_loc; 1609d6ab1dc0SHong Zhang Mat_SeqAIJ *a_loc; 1610d6ab1dc0SHong Zhang 1611d6ab1dc0SHong Zhang PetscFunctionBegin; 1612ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr); 1613d6ab1dc0SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1614d6ab1dc0SHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 1615d6ab1dc0SHong Zhang 16163cdca5ebSHong Zhang ptap = c->ap; 16178a9b8717SHong Zhang if (!ptap) SETERRQ(comm,PETSC_ERR_ARG_WRONGSTATE,"PtA cannot be reused. Do not call MatFreeIntermediateDataStructures() or use '-mat_freeintermediatedatastructures'"); 1618d6ab1dc0SHong Zhang merge = ptap->merge; 1619d6ab1dc0SHong Zhang 1620e745005bSHong Zhang /* 2) compute numeric C_seq = P_loc^T*A_loc */ 1621e745005bSHong Zhang /*------------------------------------------*/ 1622d6ab1dc0SHong Zhang /* get data from symbolic products */ 1623d6ab1dc0SHong Zhang coi = merge->coi; coj = merge->coj; 1624854ce69bSBarry Smith ierr = PetscCalloc1(coi[pon]+1,&coa);CHKERRQ(ierr); 1625d6ab1dc0SHong Zhang bi = merge->bi; bj = merge->bj; 1626d6ab1dc0SHong Zhang owners = merge->rowmap->range; 16271795a4d1SJed Brown ierr = PetscCalloc1(bi[cm]+1,&ba);CHKERRQ(ierr); 1628d6ab1dc0SHong Zhang 1629d6ab1dc0SHong Zhang /* get A_loc by taking all local rows of A */ 1630d6ab1dc0SHong Zhang A_loc = ptap->A_loc; 1631d6ab1dc0SHong Zhang ierr = MatMPIAIJGetLocalMat(A,MAT_REUSE_MATRIX,&A_loc);CHKERRQ(ierr); 1632d6ab1dc0SHong Zhang a_loc = (Mat_SeqAIJ*)(A_loc)->data; 1633d6ab1dc0SHong Zhang ai = a_loc->i; 1634d6ab1dc0SHong Zhang aj = a_loc->j; 1635d6ab1dc0SHong Zhang 1636d6ab1dc0SHong Zhang for (i=0; i<am; i++) { 1637d6ab1dc0SHong Zhang anz = ai[i+1] - ai[i]; 1638d6ab1dc0SHong Zhang adj = aj + ai[i]; 1639d6ab1dc0SHong Zhang ada = a_loc->a + ai[i]; 1640d6ab1dc0SHong Zhang 1641d6ab1dc0SHong Zhang /* 2-b) Compute Cseq = P_loc[i,:]^T*A[i,:] using outer product */ 1642e745005bSHong Zhang /*-------------------------------------------------------------*/ 1643d6ab1dc0SHong Zhang /* put the value into Co=(p->B)^T*A (off-diagonal part, send to others) */ 1644d6ab1dc0SHong Zhang pnz = po->i[i+1] - po->i[i]; 1645d6ab1dc0SHong Zhang poJ = po->j + po->i[i]; 1646d6ab1dc0SHong Zhang pA = po->a + po->i[i]; 1647d6ab1dc0SHong Zhang for (j=0; j<pnz; j++) { 1648d6ab1dc0SHong Zhang row = poJ[j]; 1649d6ab1dc0SHong Zhang cj = coj + coi[row]; 1650d6ab1dc0SHong Zhang ca = coa + coi[row]; 1651e745005bSHong Zhang /* perform sparse axpy */ 1652e745005bSHong Zhang nexta = 0; 1653d6ab1dc0SHong Zhang valtmp = pA[j]; 1654e745005bSHong Zhang for (k=0; nexta<anz; k++) { 1655e745005bSHong Zhang if (cj[k] == adj[nexta]) { 1656e745005bSHong Zhang ca[k] += valtmp*ada[nexta]; 1657e745005bSHong Zhang nexta++; 1658d6ab1dc0SHong Zhang } 1659e745005bSHong Zhang } 1660e745005bSHong Zhang ierr = PetscLogFlops(2.0*anz);CHKERRQ(ierr); 1661d6ab1dc0SHong Zhang } 1662d6ab1dc0SHong Zhang 1663d6ab1dc0SHong Zhang /* put the value into Cd (diagonal part) */ 1664d6ab1dc0SHong Zhang pnz = pd->i[i+1] - pd->i[i]; 1665d6ab1dc0SHong Zhang pdJ = pd->j + pd->i[i]; 1666d6ab1dc0SHong Zhang pA = pd->a + pd->i[i]; 1667d6ab1dc0SHong Zhang for (j=0; j<pnz; j++) { 1668d6ab1dc0SHong Zhang row = pdJ[j]; 1669d6ab1dc0SHong Zhang cj = bj + bi[row]; 1670d6ab1dc0SHong Zhang ca = ba + bi[row]; 1671e745005bSHong Zhang /* perform sparse axpy */ 1672e745005bSHong Zhang nexta = 0; 1673d6ab1dc0SHong Zhang valtmp = pA[j]; 1674e745005bSHong Zhang for (k=0; nexta<anz; k++) { 1675e745005bSHong Zhang if (cj[k] == adj[nexta]) { 1676e745005bSHong Zhang ca[k] += valtmp*ada[nexta]; 1677e745005bSHong Zhang nexta++; 1678d6ab1dc0SHong Zhang } 1679e745005bSHong Zhang } 1680e745005bSHong Zhang ierr = PetscLogFlops(2.0*anz);CHKERRQ(ierr); 1681d6ab1dc0SHong Zhang } 1682d6ab1dc0SHong Zhang } 1683d6ab1dc0SHong Zhang 1684d6ab1dc0SHong Zhang /* 3) send and recv matrix values coa */ 1685d6ab1dc0SHong Zhang /*------------------------------------*/ 1686d6ab1dc0SHong Zhang buf_ri = merge->buf_ri; 1687d6ab1dc0SHong Zhang buf_rj = merge->buf_rj; 1688d6ab1dc0SHong Zhang len_s = merge->len_s; 1689d6ab1dc0SHong Zhang ierr = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr); 1690d6ab1dc0SHong Zhang ierr = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr); 1691d6ab1dc0SHong Zhang 1692dcca6d9dSJed Brown ierr = PetscMalloc2(merge->nsend+1,&s_waits,size,&status);CHKERRQ(ierr); 1693d6ab1dc0SHong Zhang for (proc=0,k=0; proc<size; proc++) { 1694d6ab1dc0SHong Zhang if (!len_s[proc]) continue; 1695d6ab1dc0SHong Zhang i = merge->owners_co[proc]; 1696d6ab1dc0SHong Zhang ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr); 1697d6ab1dc0SHong Zhang k++; 1698d6ab1dc0SHong Zhang } 1699d6ab1dc0SHong Zhang if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);} 1700d6ab1dc0SHong Zhang if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);} 1701d6ab1dc0SHong Zhang 1702d6ab1dc0SHong Zhang ierr = PetscFree2(s_waits,status);CHKERRQ(ierr); 1703d6ab1dc0SHong Zhang ierr = PetscFree(r_waits);CHKERRQ(ierr); 1704d6ab1dc0SHong Zhang ierr = PetscFree(coa);CHKERRQ(ierr); 1705d6ab1dc0SHong Zhang 1706d6ab1dc0SHong Zhang /* 4) insert local Cseq and received values into Cmpi */ 1707d6ab1dc0SHong Zhang /*----------------------------------------------------*/ 1708dcca6d9dSJed Brown ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr); 1709d6ab1dc0SHong Zhang for (k=0; k<merge->nrecv; k++) { 1710e745005bSHong Zhang buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 1711d6ab1dc0SHong Zhang nrows = *(buf_ri_k[k]); 1712d6ab1dc0SHong Zhang nextrow[k] = buf_ri_k[k]+1; /* next row number of k-th recved i-structure */ 1713d6ab1dc0SHong Zhang nextci[k] = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure */ 1714d6ab1dc0SHong Zhang } 1715d6ab1dc0SHong Zhang 1716d6ab1dc0SHong Zhang for (i=0; i<cm; i++) { 1717d6ab1dc0SHong Zhang row = owners[rank] + i; /* global row index of C_seq */ 1718d6ab1dc0SHong Zhang bj_i = bj + bi[i]; /* col indices of the i-th row of C */ 1719d6ab1dc0SHong Zhang ba_i = ba + bi[i]; 1720d6ab1dc0SHong Zhang bnz = bi[i+1] - bi[i]; 1721d6ab1dc0SHong Zhang /* add received vals into ba */ 1722d6ab1dc0SHong Zhang for (k=0; k<merge->nrecv; k++) { /* k-th received message */ 1723d6ab1dc0SHong Zhang /* i-th row */ 1724d6ab1dc0SHong Zhang if (i == *nextrow[k]) { 1725d6ab1dc0SHong Zhang cnz = *(nextci[k]+1) - *nextci[k]; 1726d6ab1dc0SHong Zhang cj = buf_rj[k] + *(nextci[k]); 1727d6ab1dc0SHong Zhang ca = abuf_r[k] + *(nextci[k]); 1728d6ab1dc0SHong Zhang nextcj = 0; 1729d6ab1dc0SHong Zhang for (j=0; nextcj<cnz; j++) { 1730d6ab1dc0SHong Zhang if (bj_i[j] == cj[nextcj]) { /* bcol == ccol */ 1731d6ab1dc0SHong Zhang ba_i[j] += ca[nextcj++]; 1732d6ab1dc0SHong Zhang } 1733d6ab1dc0SHong Zhang } 1734d6ab1dc0SHong Zhang nextrow[k]++; nextci[k]++; 1735d6ab1dc0SHong Zhang ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr); 1736d6ab1dc0SHong Zhang } 1737d6ab1dc0SHong Zhang } 1738d6ab1dc0SHong Zhang ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr); 1739d6ab1dc0SHong Zhang } 1740d6ab1dc0SHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1741d6ab1dc0SHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1742d6ab1dc0SHong Zhang 1743d6ab1dc0SHong Zhang ierr = PetscFree(ba);CHKERRQ(ierr); 1744d6ab1dc0SHong Zhang ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr); 1745d6ab1dc0SHong Zhang ierr = PetscFree(abuf_r);CHKERRQ(ierr); 1746d6ab1dc0SHong Zhang ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr); 17478a9b8717SHong Zhang 1748*7d0a89b7SHong Zhang if (ptap->freestruct) { 17498a9b8717SHong Zhang ierr = MatFreeIntermediateDataStructures(C);CHKERRQ(ierr); 17508a9b8717SHong Zhang } 1751d6ab1dc0SHong Zhang PetscFunctionReturn(0); 1752d6ab1dc0SHong Zhang } 1753d6ab1dc0SHong Zhang 17546da69ca6SHong Zhang PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(Mat P,Mat A,PetscReal fill,Mat *C) 1755d6ab1dc0SHong Zhang { 1756d6ab1dc0SHong Zhang PetscErrorCode ierr; 1757d6ab1dc0SHong Zhang Mat Cmpi,A_loc,POt,PDt; 17583cdca5ebSHong Zhang Mat_APMPI *ptap; 17590298fd71SBarry Smith PetscFreeSpaceList free_space=NULL,current_space=NULL; 176033ba5e2fSHong Zhang Mat_MPIAIJ *p=(Mat_MPIAIJ*)P->data,*a=(Mat_MPIAIJ*)A->data,*c; 1761d6ab1dc0SHong Zhang PetscInt *pdti,*pdtj,*poti,*potj,*ptJ; 1762d6ab1dc0SHong Zhang PetscInt nnz; 1763d6ab1dc0SHong Zhang PetscInt *lnk,*owners_co,*coi,*coj,i,k,pnz,row; 1764d6ab1dc0SHong Zhang PetscInt am =A->rmap->n,pn=P->cmap->n; 1765ce94432eSBarry Smith MPI_Comm comm; 1766d6ab1dc0SHong Zhang PetscMPIInt size,rank,tagi,tagj,*len_si,*len_s,*len_ri; 1767d6ab1dc0SHong Zhang PetscInt **buf_rj,**buf_ri,**buf_ri_k; 1768d6ab1dc0SHong Zhang PetscInt len,proc,*dnz,*onz,*owners; 1769d6ab1dc0SHong Zhang PetscInt nzi,*bi,*bj; 1770d6ab1dc0SHong Zhang PetscInt nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci; 1771d6ab1dc0SHong Zhang MPI_Request *swaits,*rwaits; 1772d6ab1dc0SHong Zhang MPI_Status *sstatus,rstatus; 1773d6ab1dc0SHong Zhang Mat_Merge_SeqsToMPI *merge; 1774d6ab1dc0SHong Zhang PetscInt *ai,*aj,*Jptr,anz,*prmap=p->garray,pon,nspacedouble=0,j; 1775d6ab1dc0SHong Zhang PetscReal afill =1.0,afill_tmp; 17764b38b95cSHong Zhang PetscInt rstart = P->cmap->rstart,rmax,aN=A->cmap->N,Armax; 1777d6ab1dc0SHong Zhang PetscScalar *vals; 1778d6ab1dc0SHong Zhang Mat_SeqAIJ *a_loc,*pdt,*pot; 17790acc65b4SHong Zhang PetscTable ta; 1780b92f168fSBarry Smith MatType mtype; 1781d6ab1dc0SHong Zhang 1782d6ab1dc0SHong Zhang PetscFunctionBegin; 1783ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 1784d6ab1dc0SHong Zhang /* check if matrix local sizes are compatible */ 17856c4ed002SBarry Smith if (A->rmap->rstart != P->rmap->rstart || A->rmap->rend != P->rmap->rend) 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); 1786d6ab1dc0SHong Zhang 1787d6ab1dc0SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1788d6ab1dc0SHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 1789d6ab1dc0SHong Zhang 17903cdca5ebSHong Zhang /* create struct Mat_APMPI and attached it to C later */ 1791b00a9115SJed Brown ierr = PetscNew(&ptap);CHKERRQ(ierr); 1792d6ab1dc0SHong Zhang 1793d6ab1dc0SHong Zhang /* get A_loc by taking all local rows of A */ 1794d6ab1dc0SHong Zhang ierr = MatMPIAIJGetLocalMat(A,MAT_INITIAL_MATRIX,&A_loc);CHKERRQ(ierr); 17952205254eSKarl Rupp 1796d6ab1dc0SHong Zhang ptap->A_loc = A_loc; 1797d6ab1dc0SHong Zhang a_loc = (Mat_SeqAIJ*)(A_loc)->data; 1798d6ab1dc0SHong Zhang ai = a_loc->i; 1799d6ab1dc0SHong Zhang aj = a_loc->j; 1800d6ab1dc0SHong Zhang 1801d6ab1dc0SHong Zhang /* determine symbolic Co=(p->B)^T*A - send to others */ 1802d6ab1dc0SHong Zhang /*----------------------------------------------------*/ 1803d6ab1dc0SHong Zhang ierr = MatTransposeSymbolic_SeqAIJ(p->A,&PDt);CHKERRQ(ierr); 1804d6ab1dc0SHong Zhang pdt = (Mat_SeqAIJ*)PDt->data; 1805d6ab1dc0SHong Zhang pdti = pdt->i; pdtj = pdt->j; 1806d6ab1dc0SHong Zhang 1807d6ab1dc0SHong Zhang ierr = MatTransposeSymbolic_SeqAIJ(p->B,&POt);CHKERRQ(ierr); 1808d6ab1dc0SHong Zhang pot = (Mat_SeqAIJ*)POt->data; 1809d6ab1dc0SHong Zhang poti = pot->i; potj = pot->j; 1810d6ab1dc0SHong Zhang 1811d6ab1dc0SHong Zhang /* then, compute symbolic Co = (p->B)^T*A */ 1812d6ab1dc0SHong Zhang pon = (p->B)->cmap->n; /* total num of rows to be sent to other processors 1813d6ab1dc0SHong Zhang >= (num of nonzero rows of C_seq) - pn */ 1814854ce69bSBarry Smith ierr = PetscMalloc1(pon+1,&coi);CHKERRQ(ierr); 1815d6ab1dc0SHong Zhang coi[0] = 0; 1816d6ab1dc0SHong Zhang 1817d6ab1dc0SHong Zhang /* set initial free space to be fill*(nnz(p->B) + nnz(A)) */ 1818f91af8c7SBarry Smith nnz = PetscRealIntMultTruncate(fill,PetscIntSumTruncate(poti[pon],ai[am])); 1819a2ea699eSBarry Smith ierr = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr); 1820d6ab1dc0SHong Zhang current_space = free_space; 182119f0e57aSHong Zhang 1822d6ab1dc0SHong Zhang /* create and initialize a linked list */ 182333ba5e2fSHong Zhang ierr = PetscTableCreate(A->cmap->n + a->B->cmap->N,aN,&ta);CHKERRQ(ierr); 18244b38b95cSHong Zhang MatRowMergeMax_SeqAIJ(a_loc,am,ta); 18250acc65b4SHong Zhang ierr = PetscTableGetCount(ta,&Armax);CHKERRQ(ierr); 182633ba5e2fSHong Zhang 18270acc65b4SHong Zhang ierr = PetscLLCondensedCreate_Scalable(Armax,&lnk);CHKERRQ(ierr); 1828d6ab1dc0SHong Zhang 1829d6ab1dc0SHong Zhang for (i=0; i<pon; i++) { 1830d6ab1dc0SHong Zhang pnz = poti[i+1] - poti[i]; 1831d6ab1dc0SHong Zhang ptJ = potj + poti[i]; 1832d6ab1dc0SHong Zhang for (j=0; j<pnz; j++) { 1833d6ab1dc0SHong Zhang row = ptJ[j]; /* row of A_loc == col of Pot */ 1834d6ab1dc0SHong Zhang anz = ai[row+1] - ai[row]; 1835d6ab1dc0SHong Zhang Jptr = aj + ai[row]; 1836d6ab1dc0SHong Zhang /* add non-zero cols of AP into the sorted linked list lnk */ 1837d6ab1dc0SHong Zhang ierr = PetscLLCondensedAddSorted_Scalable(anz,Jptr,lnk);CHKERRQ(ierr); 1838d6ab1dc0SHong Zhang } 1839d6ab1dc0SHong Zhang nnz = lnk[0]; 1840d6ab1dc0SHong Zhang 1841d6ab1dc0SHong Zhang /* If free space is not available, double the total space in the list */ 1842d6ab1dc0SHong Zhang if (current_space->local_remaining<nnz) { 1843f91af8c7SBarry Smith ierr = PetscFreeSpaceGet(PetscIntSumTruncate(nnz,current_space->total_array_size),¤t_space);CHKERRQ(ierr); 1844d6ab1dc0SHong Zhang nspacedouble++; 1845d6ab1dc0SHong Zhang } 1846d6ab1dc0SHong Zhang 1847d6ab1dc0SHong Zhang /* Copy data into free space, and zero out denserows */ 1848d6ab1dc0SHong Zhang ierr = PetscLLCondensedClean_Scalable(nnz,current_space->array,lnk);CHKERRQ(ierr); 18492205254eSKarl Rupp 1850d6ab1dc0SHong Zhang current_space->array += nnz; 1851d6ab1dc0SHong Zhang current_space->local_used += nnz; 1852d6ab1dc0SHong Zhang current_space->local_remaining -= nnz; 18532205254eSKarl Rupp 1854d6ab1dc0SHong Zhang coi[i+1] = coi[i] + nnz; 1855d6ab1dc0SHong Zhang } 1856d6ab1dc0SHong Zhang 1857854ce69bSBarry Smith ierr = PetscMalloc1(coi[pon]+1,&coj);CHKERRQ(ierr); 1858d6ab1dc0SHong Zhang ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr); 18590acc65b4SHong Zhang ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr); /* must destroy to get a new one for C */ 18602205254eSKarl Rupp 1861118727c9SMark F. Adams afill_tmp = (PetscReal)coi[pon]/(poti[pon] + ai[am]+1); 1862d6ab1dc0SHong Zhang if (afill_tmp > afill) afill = afill_tmp; 1863d6ab1dc0SHong Zhang 1864d6ab1dc0SHong Zhang /* send j-array (coj) of Co to other processors */ 1865d6ab1dc0SHong Zhang /*----------------------------------------------*/ 1866d6ab1dc0SHong Zhang /* determine row ownership */ 1867b00a9115SJed Brown ierr = PetscNew(&merge);CHKERRQ(ierr); 1868d6ab1dc0SHong Zhang ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr); 18692205254eSKarl Rupp 1870d6ab1dc0SHong Zhang merge->rowmap->n = pn; 1871d6ab1dc0SHong Zhang merge->rowmap->bs = 1; 18722205254eSKarl Rupp 1873d6ab1dc0SHong Zhang ierr = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr); 1874d6ab1dc0SHong Zhang owners = merge->rowmap->range; 1875d6ab1dc0SHong Zhang 1876d6ab1dc0SHong Zhang /* determine the number of messages to send, their lengths */ 18771795a4d1SJed Brown ierr = PetscCalloc1(size,&len_si);CHKERRQ(ierr); 1878785e854fSJed Brown ierr = PetscMalloc1(size,&merge->len_s);CHKERRQ(ierr); 18792205254eSKarl Rupp 1880d6ab1dc0SHong Zhang len_s = merge->len_s; 1881d6ab1dc0SHong Zhang merge->nsend = 0; 1882d6ab1dc0SHong Zhang 1883854ce69bSBarry Smith ierr = PetscMalloc1(size+2,&owners_co);CHKERRQ(ierr); 1884d6ab1dc0SHong Zhang ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 1885d6ab1dc0SHong Zhang 1886d6ab1dc0SHong Zhang proc = 0; 1887d6ab1dc0SHong Zhang for (i=0; i<pon; i++) { 1888d6ab1dc0SHong Zhang while (prmap[i] >= owners[proc+1]) proc++; 1889d6ab1dc0SHong Zhang len_si[proc]++; /* num of rows in Co to be sent to [proc] */ 1890d6ab1dc0SHong Zhang len_s[proc] += coi[i+1] - coi[i]; 1891d6ab1dc0SHong Zhang } 1892d6ab1dc0SHong Zhang 1893d6ab1dc0SHong Zhang len = 0; /* max length of buf_si[] */ 1894d6ab1dc0SHong Zhang owners_co[0] = 0; 1895d6ab1dc0SHong Zhang for (proc=0; proc<size; proc++) { 1896d6ab1dc0SHong Zhang owners_co[proc+1] = owners_co[proc] + len_si[proc]; 1897d6ab1dc0SHong Zhang if (len_si[proc]) { 1898d6ab1dc0SHong Zhang merge->nsend++; 1899d6ab1dc0SHong Zhang len_si[proc] = 2*(len_si[proc] + 1); 1900d6ab1dc0SHong Zhang len += len_si[proc]; 1901d6ab1dc0SHong Zhang } 1902d6ab1dc0SHong Zhang } 1903d6ab1dc0SHong Zhang 1904d6ab1dc0SHong Zhang /* determine the number and length of messages to receive for coi and coj */ 19050298fd71SBarry Smith ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&merge->nrecv);CHKERRQ(ierr); 1906d6ab1dc0SHong Zhang ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr); 1907d6ab1dc0SHong Zhang 1908d6ab1dc0SHong Zhang /* post the Irecv and Isend of coj */ 1909d6ab1dc0SHong Zhang ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr); 1910d6ab1dc0SHong Zhang ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr); 1911854ce69bSBarry Smith ierr = PetscMalloc1(merge->nsend+1,&swaits);CHKERRQ(ierr); 1912d6ab1dc0SHong Zhang for (proc=0, k=0; proc<size; proc++) { 1913d6ab1dc0SHong Zhang if (!len_s[proc]) continue; 1914d6ab1dc0SHong Zhang i = owners_co[proc]; 1915d6ab1dc0SHong Zhang ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr); 1916d6ab1dc0SHong Zhang k++; 1917d6ab1dc0SHong Zhang } 1918d6ab1dc0SHong Zhang 1919d6ab1dc0SHong Zhang /* receives and sends of coj are complete */ 1920785e854fSJed Brown ierr = PetscMalloc1(size,&sstatus);CHKERRQ(ierr); 1921d6ab1dc0SHong Zhang for (i=0; i<merge->nrecv; i++) { 1922c280ed6aSJed Brown PetscMPIInt icompleted; 1923c280ed6aSJed Brown ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 1924d6ab1dc0SHong Zhang } 1925d6ab1dc0SHong Zhang ierr = PetscFree(rwaits);CHKERRQ(ierr); 1926d6ab1dc0SHong Zhang if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);} 1927d6ab1dc0SHong Zhang 19280acc65b4SHong Zhang /* add received column indices into table to update Armax */ 192933ba5e2fSHong Zhang /* Armax can be as large as aN if a P[row,:] is dense, see src/ksp/ksp/examples/tutorials/ex56.c! */ 19300acc65b4SHong Zhang for (k=0; k<merge->nrecv; k++) {/* k-th received message */ 19310acc65b4SHong Zhang Jptr = buf_rj[k]; 19320acc65b4SHong Zhang for (j=0; j<merge->len_r[k]; j++) { 19330acc65b4SHong Zhang ierr = PetscTableAdd(ta,*(Jptr+j)+1,1,INSERT_VALUES);CHKERRQ(ierr); 19340acc65b4SHong Zhang } 19350acc65b4SHong Zhang } 19360acc65b4SHong Zhang ierr = PetscTableGetCount(ta,&Armax);CHKERRQ(ierr); 193733ba5e2fSHong Zhang /* printf("Armax %d, an %d + Bn %d = %d, aN %d\n",Armax,A->cmap->n,a->B->cmap->N,A->cmap->n+a->B->cmap->N,aN); */ 19380acc65b4SHong Zhang 1939d6ab1dc0SHong Zhang /* send and recv coi */ 1940d6ab1dc0SHong Zhang /*-------------------*/ 1941d6ab1dc0SHong Zhang ierr = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr); 1942d6ab1dc0SHong Zhang ierr = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr); 1943854ce69bSBarry Smith ierr = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr); 1944d6ab1dc0SHong Zhang buf_si = buf_s; /* points to the beginning of k-th msg to be sent */ 1945d6ab1dc0SHong Zhang for (proc=0,k=0; proc<size; proc++) { 1946d6ab1dc0SHong Zhang if (!len_s[proc]) continue; 1947d6ab1dc0SHong Zhang /* form outgoing message for i-structure: 1948d6ab1dc0SHong Zhang buf_si[0]: nrows to be sent 1949d6ab1dc0SHong Zhang [1:nrows]: row index (global) 1950d6ab1dc0SHong Zhang [nrows+1:2*nrows+1]: i-structure index 1951d6ab1dc0SHong Zhang */ 1952d6ab1dc0SHong Zhang /*-------------------------------------------*/ 1953d6ab1dc0SHong Zhang nrows = len_si[proc]/2 - 1; 1954d6ab1dc0SHong Zhang buf_si_i = buf_si + nrows+1; 1955d6ab1dc0SHong Zhang buf_si[0] = nrows; 1956d6ab1dc0SHong Zhang buf_si_i[0] = 0; 1957d6ab1dc0SHong Zhang nrows = 0; 1958d6ab1dc0SHong Zhang for (i=owners_co[proc]; i<owners_co[proc+1]; i++) { 1959d6ab1dc0SHong Zhang nzi = coi[i+1] - coi[i]; 1960d6ab1dc0SHong Zhang buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */ 1961d6ab1dc0SHong Zhang buf_si[nrows+1] = prmap[i] -owners[proc]; /* local row index */ 1962d6ab1dc0SHong Zhang nrows++; 1963d6ab1dc0SHong Zhang } 1964d6ab1dc0SHong Zhang ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr); 1965d6ab1dc0SHong Zhang k++; 1966d6ab1dc0SHong Zhang buf_si += len_si[proc]; 1967d6ab1dc0SHong Zhang } 1968d6ab1dc0SHong Zhang i = merge->nrecv; 1969d6ab1dc0SHong Zhang while (i--) { 1970c280ed6aSJed Brown PetscMPIInt icompleted; 1971c280ed6aSJed Brown ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 1972d6ab1dc0SHong Zhang } 1973d6ab1dc0SHong Zhang ierr = PetscFree(rwaits);CHKERRQ(ierr); 1974d6ab1dc0SHong Zhang if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);} 1975d6ab1dc0SHong Zhang ierr = PetscFree(len_si);CHKERRQ(ierr); 1976d6ab1dc0SHong Zhang ierr = PetscFree(len_ri);CHKERRQ(ierr); 1977d6ab1dc0SHong Zhang ierr = PetscFree(swaits);CHKERRQ(ierr); 1978d6ab1dc0SHong Zhang ierr = PetscFree(sstatus);CHKERRQ(ierr); 1979d6ab1dc0SHong Zhang ierr = PetscFree(buf_s);CHKERRQ(ierr); 1980d6ab1dc0SHong Zhang 1981d6ab1dc0SHong Zhang /* compute the local portion of C (mpi mat) */ 1982d6ab1dc0SHong Zhang /*------------------------------------------*/ 1983d6ab1dc0SHong Zhang /* allocate bi array and free space for accumulating nonzero column info */ 1984854ce69bSBarry Smith ierr = PetscMalloc1(pn+1,&bi);CHKERRQ(ierr); 1985d6ab1dc0SHong Zhang bi[0] = 0; 1986d6ab1dc0SHong Zhang 1987d6ab1dc0SHong Zhang /* set initial free space to be fill*(nnz(P) + nnz(AP)) */ 1988f91af8c7SBarry Smith nnz = PetscRealIntMultTruncate(fill,PetscIntSumTruncate(pdti[pn],PetscIntSumTruncate(poti[pon],ai[am]))); 1989a2ea699eSBarry Smith ierr = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr); 1990d6ab1dc0SHong Zhang current_space = free_space; 1991d6ab1dc0SHong Zhang 1992dcca6d9dSJed Brown ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr); 1993d6ab1dc0SHong Zhang for (k=0; k<merge->nrecv; k++) { 1994d6ab1dc0SHong Zhang buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 1995d6ab1dc0SHong Zhang nrows = *buf_ri_k[k]; 1996d6ab1dc0SHong Zhang nextrow[k] = buf_ri_k[k] + 1; /* next row number of k-th recved i-structure */ 19972205254eSKarl Rupp nextci[k] = buf_ri_k[k] + (nrows + 1); /* points to the next i-structure of k-th recieved i-structure */ 1998d6ab1dc0SHong Zhang } 1999d6ab1dc0SHong Zhang 20000acc65b4SHong Zhang ierr = PetscLLCondensedCreate_Scalable(Armax,&lnk);CHKERRQ(ierr); 2001d6ab1dc0SHong Zhang ierr = MatPreallocateInitialize(comm,pn,A->cmap->n,dnz,onz);CHKERRQ(ierr); 2002d6ab1dc0SHong Zhang rmax = 0; 2003d6ab1dc0SHong Zhang for (i=0; i<pn; i++) { 2004d6ab1dc0SHong Zhang /* add pdt[i,:]*AP into lnk */ 2005d6ab1dc0SHong Zhang pnz = pdti[i+1] - pdti[i]; 2006d6ab1dc0SHong Zhang ptJ = pdtj + pdti[i]; 2007d6ab1dc0SHong Zhang for (j=0; j<pnz; j++) { 2008d6ab1dc0SHong Zhang row = ptJ[j]; /* row of AP == col of Pt */ 2009d6ab1dc0SHong Zhang anz = ai[row+1] - ai[row]; 2010d6ab1dc0SHong Zhang Jptr = aj + ai[row]; 2011d6ab1dc0SHong Zhang /* add non-zero cols of AP into the sorted linked list lnk */ 2012d6ab1dc0SHong Zhang ierr = PetscLLCondensedAddSorted_Scalable(anz,Jptr,lnk);CHKERRQ(ierr); 2013d6ab1dc0SHong Zhang } 2014d6ab1dc0SHong Zhang 2015d6ab1dc0SHong Zhang /* add received col data into lnk */ 2016d6ab1dc0SHong Zhang for (k=0; k<merge->nrecv; k++) { /* k-th received message */ 2017d6ab1dc0SHong Zhang if (i == *nextrow[k]) { /* i-th row */ 2018d6ab1dc0SHong Zhang nzi = *(nextci[k]+1) - *nextci[k]; 2019d6ab1dc0SHong Zhang Jptr = buf_rj[k] + *nextci[k]; 2020d6ab1dc0SHong Zhang ierr = PetscLLCondensedAddSorted_Scalable(nzi,Jptr,lnk);CHKERRQ(ierr); 2021d6ab1dc0SHong Zhang nextrow[k]++; nextci[k]++; 2022d6ab1dc0SHong Zhang } 2023d6ab1dc0SHong Zhang } 2024d6ab1dc0SHong Zhang nnz = lnk[0]; 2025d6ab1dc0SHong Zhang 2026d6ab1dc0SHong Zhang /* if free space is not available, make more free space */ 2027d6ab1dc0SHong Zhang if (current_space->local_remaining<nnz) { 2028f91af8c7SBarry Smith ierr = PetscFreeSpaceGet(PetscIntSumTruncate(nnz,current_space->total_array_size),¤t_space);CHKERRQ(ierr); 2029d6ab1dc0SHong Zhang nspacedouble++; 2030d6ab1dc0SHong Zhang } 2031d6ab1dc0SHong Zhang /* copy data into free space, then initialize lnk */ 2032d6ab1dc0SHong Zhang ierr = PetscLLCondensedClean_Scalable(nnz,current_space->array,lnk);CHKERRQ(ierr); 2033d6ab1dc0SHong Zhang ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr); 20342205254eSKarl Rupp 2035d6ab1dc0SHong Zhang current_space->array += nnz; 2036d6ab1dc0SHong Zhang current_space->local_used += nnz; 2037d6ab1dc0SHong Zhang current_space->local_remaining -= nnz; 20382205254eSKarl Rupp 2039d6ab1dc0SHong Zhang bi[i+1] = bi[i] + nnz; 2040d6ab1dc0SHong Zhang if (nnz > rmax) rmax = nnz; 2041d6ab1dc0SHong Zhang } 2042f671be37SHong Zhang ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr); 2043d6ab1dc0SHong Zhang 2044854ce69bSBarry Smith ierr = PetscMalloc1(bi[pn]+1,&bj);CHKERRQ(ierr); 2045d6ab1dc0SHong Zhang ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr); 2046118727c9SMark F. Adams afill_tmp = (PetscReal)bi[pn]/(pdti[pn] + poti[pon] + ai[am]+1); 2047d6ab1dc0SHong Zhang if (afill_tmp > afill) afill = afill_tmp; 2048d6ab1dc0SHong Zhang ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr); 20490acc65b4SHong Zhang ierr = PetscTableDestroy(&ta);CHKERRQ(ierr); 20500acc65b4SHong Zhang 2051d6ab1dc0SHong Zhang ierr = MatDestroy(&POt);CHKERRQ(ierr); 2052d6ab1dc0SHong Zhang ierr = MatDestroy(&PDt);CHKERRQ(ierr); 2053d6ab1dc0SHong Zhang 2054d6ab1dc0SHong Zhang /* create symbolic parallel matrix Cmpi - why cannot be assembled in Numeric part */ 2055d6ab1dc0SHong Zhang /*----------------------------------------------------------------------------------*/ 20561795a4d1SJed Brown ierr = PetscCalloc1(rmax+1,&vals);CHKERRQ(ierr); 2057d6ab1dc0SHong Zhang 2058d6ab1dc0SHong Zhang ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr); 2059d6ab1dc0SHong Zhang ierr = MatSetSizes(Cmpi,pn,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 206033d57670SJed Brown ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(A->cmap->bs));CHKERRQ(ierr); 2061b92f168fSBarry Smith ierr = MatGetType(A,&mtype);CHKERRQ(ierr); 2062b92f168fSBarry Smith ierr = MatSetType(Cmpi,mtype);CHKERRQ(ierr); 2063d6ab1dc0SHong Zhang ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr); 2064d6ab1dc0SHong Zhang ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 2065d6ab1dc0SHong Zhang ierr = MatSetBlockSize(Cmpi,1);CHKERRQ(ierr); 2066d6ab1dc0SHong Zhang for (i=0; i<pn; i++) { 2067d6ab1dc0SHong Zhang row = i + rstart; 2068d6ab1dc0SHong Zhang nnz = bi[i+1] - bi[i]; 2069d6ab1dc0SHong Zhang Jptr = bj + bi[i]; 2070d6ab1dc0SHong Zhang ierr = MatSetValues(Cmpi,1,&row,nnz,Jptr,vals,INSERT_VALUES);CHKERRQ(ierr); 2071d6ab1dc0SHong Zhang } 2072d6ab1dc0SHong Zhang ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2073d6ab1dc0SHong Zhang ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2074d6ab1dc0SHong Zhang ierr = PetscFree(vals);CHKERRQ(ierr); 2075d6ab1dc0SHong Zhang 2076d6ab1dc0SHong Zhang merge->bi = bi; 2077d6ab1dc0SHong Zhang merge->bj = bj; 2078d6ab1dc0SHong Zhang merge->coi = coi; 2079d6ab1dc0SHong Zhang merge->coj = coj; 2080d6ab1dc0SHong Zhang merge->buf_ri = buf_ri; 2081d6ab1dc0SHong Zhang merge->buf_rj = buf_rj; 2082d6ab1dc0SHong Zhang merge->owners_co = owners_co; 2083d6ab1dc0SHong Zhang 2084d6ab1dc0SHong Zhang /* attach the supporting struct to Cmpi for reuse */ 2085d6ab1dc0SHong Zhang c = (Mat_MPIAIJ*)Cmpi->data; 20862205254eSKarl Rupp 20873cdca5ebSHong Zhang c->ap = ptap; 20880298fd71SBarry Smith ptap->api = NULL; 20890298fd71SBarry Smith ptap->apj = NULL; 2090d6ab1dc0SHong Zhang ptap->merge = merge; 20910298fd71SBarry Smith ptap->apa = NULL; 2092a560ef98SHong Zhang ptap->destroy = Cmpi->ops->destroy; 2093a560ef98SHong Zhang ptap->duplicate = Cmpi->ops->duplicate; 2094a560ef98SHong Zhang 2095a560ef98SHong Zhang Cmpi->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ; 2096a560ef98SHong Zhang Cmpi->ops->destroy = MatDestroy_MPIAIJ_PtAP; 20974624976aSHong Zhang Cmpi->ops->freeintermediatedatastructures = MatFreeIntermediateDataStructures_MPIAIJ_AP; 2098f96d37fcSHong Zhang 2099f96d37fcSHong Zhang *C = Cmpi; 2100f96d37fcSHong Zhang #if defined(PETSC_USE_INFO) 2101f96d37fcSHong Zhang if (bi[pn] != 0) { 210257622a8eSBarry Smith ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr); 210357622a8eSBarry Smith ierr = PetscInfo1(Cmpi,"Use MatTransposeMatMult(A,B,MatReuse,%g,&C) for best performance.\n",(double)afill);CHKERRQ(ierr); 2104f96d37fcSHong Zhang } else { 2105f96d37fcSHong Zhang ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr); 2106f96d37fcSHong Zhang } 2107f96d37fcSHong Zhang #endif 2108f96d37fcSHong Zhang PetscFunctionReturn(0); 2109f96d37fcSHong Zhang } 2110