182412ba7SHong Zhang 242c57489SHong Zhang /* 342c57489SHong Zhang Defines projective product routines where A is a MPIAIJ matrix 442c57489SHong Zhang C = P^T * A * P 542c57489SHong Zhang */ 642c57489SHong Zhang 7c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/aij.h> /*I "petscmat.h" I*/ 8c6db04a5SJed Brown #include <../src/mat/utils/freespace.h> 9c6db04a5SJed Brown #include <../src/mat/impls/aij/mpi/mpiaij.h> 10c6db04a5SJed Brown #include <petscbt.h> 118563dfccSBarry Smith #include <petsctime.h> 1242c57489SHong Zhang 139516a85cSHong Zhang #define PTAP_PROFILE 1424ecddacSHong Zhang 1509573ac7SBarry Smith extern PetscErrorCode MatDestroy_MPIAIJ(Mat); 16cc6cb787SHong Zhang #undef __FUNCT__ 17f8487c73SHong Zhang #define __FUNCT__ "MatDestroy_MPIAIJ_PtAP" 18f8487c73SHong Zhang PetscErrorCode MatDestroy_MPIAIJ_PtAP(Mat A) 19cc6cb787SHong Zhang { 20cc6cb787SHong Zhang PetscErrorCode ierr; 21f8487c73SHong Zhang Mat_MPIAIJ *a =(Mat_MPIAIJ*)A->data; 22f8487c73SHong Zhang Mat_PtAPMPI *ptap=a->ptap; 23cc6cb787SHong Zhang 24cc6cb787SHong Zhang PetscFunctionBegin; 25f8487c73SHong Zhang if (ptap) { 26c8b0795cSMark F. Adams Mat_Merge_SeqsToMPI *merge=ptap->merge; 27b7f45c76SHong Zhang ierr = PetscFree2(ptap->startsj_s,ptap->startsj_r);CHKERRQ(ierr); 28f8487c73SHong Zhang ierr = PetscFree(ptap->bufa);CHKERRQ(ierr); 29a1a4e84aSHong Zhang ierr = MatDestroy(&ptap->P_loc);CHKERRQ(ierr); 30a1a4e84aSHong Zhang ierr = MatDestroy(&ptap->P_oth);CHKERRQ(ierr); 31c5af039cSHong Zhang ierr = MatDestroy(&ptap->A_loc);CHKERRQ(ierr); /* used by MatTransposeMatMult() */ 320d3441aeSHong Zhang 3305d62848SHong Zhang ierr = MatDestroy(&ptap->Rd);CHKERRQ(ierr); 3405d62848SHong Zhang ierr = MatDestroy(&ptap->Ro);CHKERRQ(ierr); 35445158ffSHong Zhang 36445158ffSHong Zhang ierr = PetscFree(ptap->api);CHKERRQ(ierr); 37445158ffSHong Zhang ierr = PetscFree2(ptap->apj,ptap->apv);CHKERRQ(ierr); 380d3441aeSHong Zhang ierr = MatDestroy(&ptap->AP_loc);CHKERRQ(ierr); 392259aa2eSHong Zhang ierr = MatDestroy(&ptap->C_loc);CHKERRQ(ierr); 402259aa2eSHong Zhang ierr = MatDestroy(&ptap->C_oth);CHKERRQ(ierr); 410d3441aeSHong Zhang 42c5af039cSHong Zhang if (ptap->api) {ierr = PetscFree(ptap->api);CHKERRQ(ierr);} 43c5af039cSHong Zhang if (ptap->apj) {ierr = PetscFree(ptap->apj);CHKERRQ(ierr);} 44414894bdSHong Zhang if (ptap->apa) {ierr = PetscFree(ptap->apa);CHKERRQ(ierr);} 45c8b0795cSMark F. Adams if (merge) { 46cc6cb787SHong Zhang ierr = PetscFree(merge->id_r);CHKERRQ(ierr); 47cc6cb787SHong Zhang ierr = PetscFree(merge->len_s);CHKERRQ(ierr); 48cc6cb787SHong Zhang ierr = PetscFree(merge->len_r);CHKERRQ(ierr); 49cc6cb787SHong Zhang ierr = PetscFree(merge->bi);CHKERRQ(ierr); 50cc6cb787SHong Zhang ierr = PetscFree(merge->bj);CHKERRQ(ierr); 51c05d87d6SBarry Smith ierr = PetscFree(merge->buf_ri[0]);CHKERRQ(ierr); 52cc6cb787SHong Zhang ierr = PetscFree(merge->buf_ri);CHKERRQ(ierr); 53c05d87d6SBarry Smith ierr = PetscFree(merge->buf_rj[0]);CHKERRQ(ierr); 54cc6cb787SHong Zhang ierr = PetscFree(merge->buf_rj);CHKERRQ(ierr); 55445158ffSHong Zhang ierr = PetscFree(merge->coi);CHKERRQ(ierr); 56445158ffSHong Zhang ierr = PetscFree(merge->coj);CHKERRQ(ierr); 5705b42c5fSBarry Smith ierr = PetscFree(merge->owners_co);CHKERRQ(ierr); 586bf464f9SBarry Smith ierr = PetscLayoutDestroy(&merge->rowmap);CHKERRQ(ierr); 59dce485f0SBarry Smith ierr = merge->destroy(A);CHKERRQ(ierr); 60f8487c73SHong Zhang ierr = PetscFree(ptap->merge);CHKERRQ(ierr); 61445158ffSHong Zhang } else { 62445158ffSHong Zhang ierr = MatDestroy_MPIAIJ(A);CHKERRQ(ierr); 63bf0cc555SLisandro Dalcin } 64f8487c73SHong Zhang ierr = PetscFree(ptap);CHKERRQ(ierr); 65c8b0795cSMark F. Adams } 66cc6cb787SHong Zhang PetscFunctionReturn(0); 67cc6cb787SHong Zhang } 68cc6cb787SHong Zhang 69cc6cb787SHong Zhang #undef __FUNCT__ 70cc6cb787SHong Zhang #define __FUNCT__ "MatDuplicate_MPIAIJ_MatPtAP" 71f0c0a3a6SBarry Smith PetscErrorCode MatDuplicate_MPIAIJ_MatPtAP(Mat A, MatDuplicateOption op, Mat *M) 72f0c0a3a6SBarry Smith { 73cc6cb787SHong Zhang PetscErrorCode ierr; 7411a6856fSHong Zhang Mat_MPIAIJ *a = (Mat_MPIAIJ*)A->data; 7511a6856fSHong Zhang Mat_PtAPMPI *ptap = a->ptap; 7611a6856fSHong Zhang Mat_Merge_SeqsToMPI *merge = ptap->merge; 77cc6cb787SHong Zhang 78cc6cb787SHong Zhang PetscFunctionBegin; 79dce485f0SBarry Smith ierr = (*merge->duplicate)(A,op,M);CHKERRQ(ierr); 802205254eSKarl Rupp 8111a6856fSHong Zhang (*M)->ops->destroy = merge->destroy; 8211a6856fSHong Zhang (*M)->ops->duplicate = merge->duplicate; 83cc6cb787SHong Zhang PetscFunctionReturn(0); 84cc6cb787SHong Zhang } 85cc6cb787SHong Zhang 8642c57489SHong Zhang #undef __FUNCT__ 87cf3ca8ceSHong Zhang #define __FUNCT__ "MatPtAP_MPIAIJ_MPIAIJ" 88cf3ca8ceSHong Zhang PetscErrorCode MatPtAP_MPIAIJ_MPIAIJ(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C) 8942c57489SHong Zhang { 9042c57489SHong Zhang PetscErrorCode ierr; 910d3441aeSHong Zhang PetscBool newalg=PETSC_FALSE; 92*67a12041SHong Zhang MPI_Comm comm; 9342c57489SHong Zhang 9442c57489SHong Zhang PetscFunctionBegin; 95*67a12041SHong Zhang /* check if matrix local sizes are compatible */ 96*67a12041SHong Zhang ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 97*67a12041SHong Zhang if (A->rmap->rstart != P->rmap->rstart || A->rmap->rend != P->rmap->rend) { 98*67a12041SHong Zhang SETERRQ4(comm,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, Arow (%D, %D) != Prow (%D,%D)",A->rmap->rstart,A->rmap->rend,P->rmap->rstart,P->rmap->rend); 99*67a12041SHong Zhang } 100*67a12041SHong Zhang if (A->cmap->rstart != P->rmap->rstart || A->cmap->rend != P->rmap->rend) { 101*67a12041SHong Zhang SETERRQ4(comm,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, Acol (%D, %D) != Prow (%D,%D)",A->cmap->rstart,A->cmap->rend,P->rmap->rstart,P->rmap->rend); 102*67a12041SHong Zhang } 103*67a12041SHong Zhang 1040d3441aeSHong Zhang ierr = PetscOptionsGetBool(NULL,"-matptap_new",&newalg,NULL);CHKERRQ(ierr); 105cf3ca8ceSHong Zhang if (scall == MAT_INITIAL_MATRIX) { 1063ff4c91cSHong Zhang ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 1070d3441aeSHong Zhang if (newalg) { 1080d3441aeSHong Zhang ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ_new(A,P,fill,C);CHKERRQ(ierr); 1090d3441aeSHong Zhang } else { 110cf3ca8ceSHong Zhang ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ(A,P,fill,C);CHKERRQ(ierr); 1110d3441aeSHong Zhang } 1123ff4c91cSHong Zhang ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 1137a7894deSKris Buschelman } 1143ff4c91cSHong Zhang ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr); 1150d3441aeSHong Zhang if (newalg) { 1160d3441aeSHong Zhang ierr = MatPtAPNumeric_MPIAIJ_MPIAIJ_new(A,P,*C);CHKERRQ(ierr); 1170d3441aeSHong Zhang } else { 118cf3ca8ceSHong Zhang ierr = MatPtAPNumeric_MPIAIJ_MPIAIJ(A,P,*C);CHKERRQ(ierr); 1190d3441aeSHong Zhang } 1203ff4c91cSHong Zhang ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr); 12142c57489SHong Zhang PetscFunctionReturn(0); 12242c57489SHong Zhang } 12342c57489SHong Zhang 12442c57489SHong Zhang #undef __FUNCT__ 1250d3441aeSHong Zhang #define __FUNCT__ "MatPtAPSymbolic_MPIAIJ_MPIAIJ_new" 1260d3441aeSHong Zhang PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ_new(Mat A,Mat P,PetscReal fill,Mat *C) 1270d3441aeSHong Zhang { 1280d3441aeSHong Zhang PetscErrorCode ierr; 1290d3441aeSHong Zhang Mat_PtAPMPI *ptap; 130*67a12041SHong Zhang Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data; 1312259aa2eSHong Zhang Mat_MPIAIJ *c; 1322259aa2eSHong Zhang MPI_Comm comm; 1332259aa2eSHong Zhang PetscMPIInt size,rank; 134445158ffSHong Zhang PetscLogDouble t0,t1,t11,t12,t2,t3,t4; 13515a3b8e2SHong Zhang Mat Cmpi; 13615a3b8e2SHong Zhang PetscFreeSpaceList free_space=NULL,current_space=NULL; 137*67a12041SHong Zhang PetscInt am=A->rmap->n,pm=P->rmap->n; 138*67a12041SHong Zhang // PetscInt *pdti,*pdtj,*ptJ,nnz; 139*67a12041SHong Zhang PetscInt *lnk,i,k,pnz,row,nnz; 140445158ffSHong Zhang PetscInt pN=P->cmap->N,pn=P->cmap->n; 14115a3b8e2SHong Zhang PetscBT lnkbt; 14215a3b8e2SHong Zhang PetscMPIInt tagi,tagj,*len_si,*len_s,*len_ri,icompleted=0; 14315a3b8e2SHong Zhang PetscInt **buf_rj,**buf_ri,**buf_ri_k; 14415a3b8e2SHong Zhang PetscInt len,proc,*dnz,*onz,*owners; 145*67a12041SHong Zhang PetscInt nzi,nspacedouble; 14615a3b8e2SHong Zhang PetscInt nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci; 14715a3b8e2SHong Zhang MPI_Request *swaits,*rwaits; 14815a3b8e2SHong Zhang MPI_Status *sstatus,rstatus; 149445158ffSHong Zhang PetscLayout rowmap; 150445158ffSHong Zhang PetscInt *owners_co,*coi,*coj; /* i and j array of (p->B)^T*A*P - used in the communication */ 151445158ffSHong Zhang PetscInt nsend,nrecv; 152445158ffSHong Zhang PetscMPIInt *len_r,*id_r; /* array of length of comm->size, store send/recv matrix values */ 153*67a12041SHong Zhang PetscInt *api,*apj,*Jptr,apnz,*prmap=p->garray,con,j,ap_rmax=0; 154*67a12041SHong Zhang PetscInt rmax,*aj,*ai,*pi; 155*67a12041SHong Zhang Mat_SeqAIJ *p_loc,*p_oth,*ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*c_loc,*c_oth; 156*67a12041SHong Zhang PetscScalar *apv; 157*67a12041SHong Zhang 158*67a12041SHong Zhang PetscFunctionBegin; 159*67a12041SHong Zhang ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 160*67a12041SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 161*67a12041SHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 16215a3b8e2SHong Zhang 16380bb4639SHong Zhang ierr = PetscTime(&t0);CHKERRQ(ierr); 16415a3b8e2SHong Zhang /* create struct Mat_PtAPMPI and attached it to C later */ 16515a3b8e2SHong Zhang ierr = PetscNew(&ptap);CHKERRQ(ierr); 16615a3b8e2SHong Zhang ptap->reuse = MAT_INITIAL_MATRIX; 16715a3b8e2SHong Zhang 16815a3b8e2SHong Zhang /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */ 16915a3b8e2SHong Zhang ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 17015a3b8e2SHong Zhang /* get P_loc by taking all local rows of P */ 17115a3b8e2SHong Zhang ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 17215a3b8e2SHong Zhang 173*67a12041SHong Zhang /* (0) compute Rd = Pd^T, Ro = Po^T */ 174*67a12041SHong Zhang /* --------------------------------- */ 175de817e96SHong Zhang ierr = MatTranspose_SeqAIJ(p->A,MAT_INITIAL_MATRIX,&ptap->Rd);CHKERRQ(ierr); 176de817e96SHong Zhang ierr = MatTranspose_SeqAIJ(p->B,MAT_INITIAL_MATRIX,&ptap->Ro);CHKERRQ(ierr); 177445158ffSHong Zhang ierr = PetscTime(&t11);CHKERRQ(ierr); 17815a3b8e2SHong Zhang 179*67a12041SHong Zhang /* (1) compute symbolic AP = A_loc*P = Ad*P_loc + Ao*P_oth (api,apj) */ 180*67a12041SHong Zhang /* ----------------------------------------------------------------- */ 181*67a12041SHong Zhang p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 182*67a12041SHong Zhang p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 183445158ffSHong Zhang 184445158ffSHong Zhang /* create and initialize a linked list */ 185445158ffSHong Zhang ierr = PetscLLCondensedCreate(pN,pN,&lnk,&lnkbt);CHKERRQ(ierr); 186445158ffSHong Zhang 187445158ffSHong Zhang /* Initial FreeSpace size is fill*(nnz(A) + nnz(P)) -OOM for ex56, np=8k on Intrepid! */ 188*67a12041SHong Zhang ierr = PetscFreeSpaceGet((PetscInt)(fill*(ad->i[am]+ao->i[am]+p_loc->i[pm])),&free_space);CHKERRQ(ierr); 189445158ffSHong Zhang current_space = free_space; 190*67a12041SHong Zhang nspacedouble = 0; 191*67a12041SHong Zhang 192*67a12041SHong Zhang ierr = PetscMalloc1(am+1,&api);CHKERRQ(ierr); 193*67a12041SHong Zhang api[0] = 0; 194445158ffSHong Zhang for (i=0; i<am; i++) { 195*67a12041SHong Zhang /* diagonal portion: Ad[i,:]*P */ 196*67a12041SHong Zhang ai = ad->i; pi = p_loc->i; 197*67a12041SHong Zhang nzi = ai[i+1] - ai[i]; 198*67a12041SHong Zhang aj = ad->j + ai[i]; 199445158ffSHong Zhang for (j=0; j<nzi; j++) { 200445158ffSHong Zhang row = aj[j]; 201*67a12041SHong Zhang pnz = pi[row+1] - pi[row]; 202*67a12041SHong Zhang Jptr = p_loc->j + pi[row]; 203445158ffSHong Zhang /* add non-zero cols of P into the sorted linked list lnk */ 204445158ffSHong Zhang ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 205445158ffSHong Zhang } 206*67a12041SHong Zhang /* off-diagonal portion: Ao[i,:]*P */ 207*67a12041SHong Zhang ai = ao->i; pi = p_oth->i; 208*67a12041SHong Zhang nzi = ai[i+1] - ai[i]; 209*67a12041SHong Zhang aj = ao->j + ai[i]; 210445158ffSHong Zhang for (j=0; j<nzi; j++) { 211445158ffSHong Zhang row = aj[j]; 212*67a12041SHong Zhang pnz = pi[row+1] - pi[row]; 213*67a12041SHong Zhang Jptr = p_oth->j + pi[row]; 214445158ffSHong Zhang ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 215445158ffSHong Zhang } 216445158ffSHong Zhang apnz = lnk[0]; 217445158ffSHong Zhang api[i+1] = api[i] + apnz; 218445158ffSHong Zhang if (ap_rmax < apnz) ap_rmax = apnz; 219445158ffSHong Zhang 220445158ffSHong Zhang /* if free space is not available, double the total space in the list */ 221445158ffSHong Zhang if (current_space->local_remaining<apnz) { 222445158ffSHong Zhang ierr = PetscFreeSpaceGet(apnz+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 223445158ffSHong Zhang nspacedouble++; 224445158ffSHong Zhang } 225445158ffSHong Zhang 226445158ffSHong Zhang /* Copy data into free space, then initialize lnk */ 227445158ffSHong Zhang ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 228445158ffSHong Zhang 229445158ffSHong Zhang current_space->array += apnz; 230445158ffSHong Zhang current_space->local_used += apnz; 231445158ffSHong Zhang current_space->local_remaining -= apnz; 232445158ffSHong Zhang } 233445158ffSHong Zhang /* Allocate space for apj, initialize apj, and */ 234445158ffSHong Zhang /* destroy list of free space and other temporary array(s) */ 235445158ffSHong Zhang ierr = PetscMalloc2(api[am],&apj,api[am],&apv);CHKERRQ(ierr); 236445158ffSHong Zhang ierr = PetscFreeSpaceContiguous(&free_space,apj);CHKERRQ(ierr); 237445158ffSHong Zhang //afill_tmp = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1); 238445158ffSHong Zhang //if (afill_tmp > afill) afill = afill_tmp; 239445158ffSHong Zhang 240445158ffSHong Zhang ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,am,pN,api,apj,apv,&ptap->AP_loc);CHKERRQ(ierr); 241445158ffSHong Zhang ierr = PetscTime(&t12);CHKERRQ(ierr); 24215a3b8e2SHong Zhang 243de817e96SHong Zhang /* (2) compute C_loc = Rd*AP_loc, Co = Ro*AP_loc */ 244*67a12041SHong Zhang /* ---------------------------------------------- */ 245*67a12041SHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Rd,ptap->AP_loc,fill,&ptap->C_loc);CHKERRQ(ierr); 246*67a12041SHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Ro,ptap->AP_loc,fill,&ptap->C_oth);CHKERRQ(ierr); 247*67a12041SHong Zhang c_loc = (Mat_SeqAIJ*)ptap->C_loc->data; 248*67a12041SHong Zhang c_oth = (Mat_SeqAIJ*)ptap->C_oth->data; 24980bb4639SHong Zhang ierr = PetscTime(&t1);CHKERRQ(ierr); 25015a3b8e2SHong Zhang 251*67a12041SHong Zhang /* (3) send coj of C_oth to other processors */ 252*67a12041SHong Zhang /* ------------------------------------------ */ 253445158ffSHong Zhang ierr = PetscCalloc1(size,&len_s);CHKERRQ(ierr); 25415a3b8e2SHong Zhang 25515a3b8e2SHong Zhang /* determine row ownership */ 256445158ffSHong Zhang ierr = PetscLayoutCreate(comm,&rowmap);CHKERRQ(ierr); 257445158ffSHong Zhang rowmap->n = pn; 258445158ffSHong Zhang rowmap->bs = 1; 25915a3b8e2SHong Zhang 260445158ffSHong Zhang ierr = PetscLayoutSetUp(rowmap);CHKERRQ(ierr); 261445158ffSHong Zhang owners = rowmap->range; 26215a3b8e2SHong Zhang 26315a3b8e2SHong Zhang /* determine the number of messages to send, their lengths */ 264*67a12041SHong Zhang ierr = PetscMalloc3(size,&len_si,size,&sstatus,size+2,&owners_co);CHKERRQ(ierr); 26515a3b8e2SHong Zhang ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 26615a3b8e2SHong Zhang 267*67a12041SHong Zhang c_oth = (Mat_SeqAIJ*)ptap->C_oth->data; 268*67a12041SHong Zhang coi = c_oth->i; coj = c_oth->j; 269*67a12041SHong Zhang con = ptap->C_oth->rmap->n; 27015a3b8e2SHong Zhang proc = 0; 271*67a12041SHong Zhang for (i=0; i<con; i++) { 27215a3b8e2SHong Zhang while (prmap[i] >= owners[proc+1]) proc++; 27315a3b8e2SHong Zhang len_si[proc]++; /* num of rows in Co(=Pt*AP) to be sent to [proc] */ 27415a3b8e2SHong Zhang len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */ 27515a3b8e2SHong Zhang } 27615a3b8e2SHong Zhang 27715a3b8e2SHong Zhang len = 0; /* max length of buf_si[], see (4) */ 27815a3b8e2SHong Zhang owners_co[0] = 0; 279*67a12041SHong Zhang nsend = 0; 28015a3b8e2SHong Zhang for (proc=0; proc<size; proc++) { 28115a3b8e2SHong Zhang owners_co[proc+1] = owners_co[proc] + len_si[proc]; 28215a3b8e2SHong Zhang if (len_s[proc]) { 283445158ffSHong Zhang nsend++; 28415a3b8e2SHong Zhang len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */ 28515a3b8e2SHong Zhang len += len_si[proc]; 28615a3b8e2SHong Zhang } 28715a3b8e2SHong Zhang } 28815a3b8e2SHong Zhang 28915a3b8e2SHong Zhang /* determine the number and length of messages to receive for coi and coj */ 290445158ffSHong Zhang ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&nrecv);CHKERRQ(ierr); 291445158ffSHong Zhang ierr = PetscGatherMessageLengths2(comm,nsend,nrecv,len_s,len_si,&id_r,&len_r,&len_ri);CHKERRQ(ierr); 29215a3b8e2SHong Zhang 29315a3b8e2SHong Zhang /* post the Irecv and Isend of coj */ 29415a3b8e2SHong Zhang ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr); 295445158ffSHong Zhang ierr = PetscPostIrecvInt(comm,tagj,nrecv,id_r,len_r,&buf_rj,&rwaits);CHKERRQ(ierr); 296445158ffSHong Zhang ierr = PetscMalloc1(nsend+1,&swaits);CHKERRQ(ierr); 29715a3b8e2SHong Zhang for (proc=0, k=0; proc<size; proc++) { 29815a3b8e2SHong Zhang if (!len_s[proc]) continue; 29915a3b8e2SHong Zhang i = owners_co[proc]; 30015a3b8e2SHong Zhang ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr); 30115a3b8e2SHong Zhang k++; 30215a3b8e2SHong Zhang } 30315a3b8e2SHong Zhang 30415a3b8e2SHong Zhang /* receives and sends of coj are complete */ 305445158ffSHong Zhang for (i=0; i<nrecv; i++) { 306445158ffSHong Zhang ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 30715a3b8e2SHong Zhang } 30815a3b8e2SHong Zhang ierr = PetscFree(rwaits);CHKERRQ(ierr); 309445158ffSHong Zhang if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);} 31015a3b8e2SHong Zhang 31115a3b8e2SHong Zhang /* (4) send and recv coi */ 31215a3b8e2SHong Zhang /*-----------------------*/ 31315a3b8e2SHong Zhang ierr = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr); 314445158ffSHong Zhang ierr = PetscPostIrecvInt(comm,tagi,nrecv,id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr); 31515a3b8e2SHong Zhang ierr = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr); 31615a3b8e2SHong Zhang buf_si = buf_s; /* points to the beginning of k-th msg to be sent */ 31715a3b8e2SHong Zhang for (proc=0,k=0; proc<size; proc++) { 31815a3b8e2SHong Zhang if (!len_s[proc]) continue; 31915a3b8e2SHong Zhang /* form outgoing message for i-structure: 32015a3b8e2SHong Zhang buf_si[0]: nrows to be sent 32115a3b8e2SHong Zhang [1:nrows]: row index (global) 32215a3b8e2SHong Zhang [nrows+1:2*nrows+1]: i-structure index 32315a3b8e2SHong Zhang */ 32415a3b8e2SHong Zhang /*-------------------------------------------*/ 32515a3b8e2SHong Zhang nrows = len_si[proc]/2 - 1; /* num of rows in Co to be sent to [proc] */ 32615a3b8e2SHong Zhang buf_si_i = buf_si + nrows+1; 32715a3b8e2SHong Zhang buf_si[0] = nrows; 32815a3b8e2SHong Zhang buf_si_i[0] = 0; 32915a3b8e2SHong Zhang nrows = 0; 33015a3b8e2SHong Zhang for (i=owners_co[proc]; i<owners_co[proc+1]; i++) { 33115a3b8e2SHong Zhang nzi = coi[i+1] - coi[i]; 33215a3b8e2SHong Zhang buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */ 33315a3b8e2SHong Zhang buf_si[nrows+1] = prmap[i] -owners[proc]; /* local row index */ 33415a3b8e2SHong Zhang nrows++; 33515a3b8e2SHong Zhang } 33615a3b8e2SHong Zhang ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr); 33715a3b8e2SHong Zhang k++; 33815a3b8e2SHong Zhang buf_si += len_si[proc]; 33915a3b8e2SHong Zhang } 340445158ffSHong Zhang i = nrecv; 34115a3b8e2SHong Zhang while (i--) { 342445158ffSHong Zhang ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 34315a3b8e2SHong Zhang } 34415a3b8e2SHong Zhang ierr = PetscFree(rwaits);CHKERRQ(ierr); 345445158ffSHong Zhang if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);} 34615a3b8e2SHong Zhang 347*67a12041SHong Zhang ierr = PetscFree3(len_si,sstatus,owners_co);CHKERRQ(ierr); 34815a3b8e2SHong Zhang ierr = PetscFree(len_ri);CHKERRQ(ierr); 34915a3b8e2SHong Zhang ierr = PetscFree(swaits);CHKERRQ(ierr); 35015a3b8e2SHong Zhang ierr = PetscFree(buf_s);CHKERRQ(ierr); 35115a3b8e2SHong Zhang 35280bb4639SHong Zhang ierr = PetscTime(&t2);CHKERRQ(ierr); 35380bb4639SHong Zhang 354*67a12041SHong Zhang /* (5) compute the local portion of Cmpi */ 355*67a12041SHong Zhang /* ------------------------------------------ */ 35615a3b8e2SHong Zhang /* set initial free space to be fill*(nnz(P) + nnz(AP)) */ 357445158ffSHong Zhang //nnz = fill*(pi_loc[pm] + api[am]); 358445158ffSHong Zhang nnz = pN; //new - nonscalable!!! 359*67a12041SHong Zhang //if (!rank) printf("pN %d, C_loc->rmax %d, C_oth->rmax %d\n",pN,c_loc->rmax,c_oth->rmax); 36015a3b8e2SHong Zhang ierr = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr); 36115a3b8e2SHong Zhang current_space = free_space; 36215a3b8e2SHong Zhang 363445158ffSHong Zhang ierr = PetscMalloc3(nrecv,&buf_ri_k,nrecv,&nextrow,nrecv,&nextci);CHKERRQ(ierr); 364445158ffSHong Zhang for (k=0; k<nrecv; k++) { 36515a3b8e2SHong Zhang buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 36615a3b8e2SHong Zhang nrows = *buf_ri_k[k]; 36715a3b8e2SHong Zhang nextrow[k] = buf_ri_k[k] + 1; /* next row number of k-th recved i-structure */ 36815a3b8e2SHong Zhang nextci[k] = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure */ 36915a3b8e2SHong Zhang } 37015a3b8e2SHong Zhang ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr); 371445158ffSHong Zhang 37215a3b8e2SHong Zhang rmax = 0; 37315a3b8e2SHong Zhang for (i=0; i<pn; i++) { 374*67a12041SHong Zhang /* add C_loc into Cmpi */ 375*67a12041SHong Zhang nzi = c_loc->i[i+1] - c_loc->i[i]; 376*67a12041SHong Zhang Jptr = c_loc->j + c_loc->i[i]; 377*67a12041SHong Zhang ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr); 37815a3b8e2SHong Zhang 37915a3b8e2SHong Zhang /* add received col data into lnk */ 380445158ffSHong Zhang for (k=0; k<nrecv; k++) { /* k-th received message */ 38115a3b8e2SHong Zhang if (i == *nextrow[k]) { /* i-th row */ 38215a3b8e2SHong Zhang nzi = *(nextci[k]+1) - *nextci[k]; 38315a3b8e2SHong Zhang Jptr = buf_rj[k] + *nextci[k]; 38415a3b8e2SHong Zhang ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr); 38515a3b8e2SHong Zhang nextrow[k]++; nextci[k]++; 38615a3b8e2SHong Zhang } 38715a3b8e2SHong Zhang } 38815a3b8e2SHong Zhang nnz = lnk[0]; 389445158ffSHong Zhang #if 0 39015a3b8e2SHong Zhang /* if free space is not available, make more free space */ 39115a3b8e2SHong Zhang if (current_space->local_remaining<nnz) { 39215a3b8e2SHong Zhang ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 39315a3b8e2SHong Zhang nspacedouble++; 39415a3b8e2SHong Zhang } 395445158ffSHong Zhang #endif 39615a3b8e2SHong Zhang /* copy data into free space, then initialize lnk */ 39715a3b8e2SHong Zhang ierr = PetscLLCondensedClean(pN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 39815a3b8e2SHong Zhang ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr); 39915a3b8e2SHong Zhang 40015a3b8e2SHong Zhang if (nnz > rmax) rmax = nnz; 40115a3b8e2SHong Zhang } 40215a3b8e2SHong Zhang ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr); 40315a3b8e2SHong Zhang ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 404445158ffSHong Zhang ierr = PetscFreeSpaceDestroy(free_space);CHKERRQ(ierr); 40515a3b8e2SHong Zhang 40680bb4639SHong Zhang ierr = PetscTime(&t3);CHKERRQ(ierr); 40780bb4639SHong Zhang 40815a3b8e2SHong Zhang /* (6) create symbolic parallel matrix Cmpi */ 40915a3b8e2SHong Zhang /*------------------------------------------*/ 41015a3b8e2SHong Zhang ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr); 41115a3b8e2SHong Zhang ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 41215a3b8e2SHong Zhang ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr); 41315a3b8e2SHong Zhang ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr); 41415a3b8e2SHong Zhang ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr); 41515a3b8e2SHong Zhang ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 41615a3b8e2SHong Zhang 417445158ffSHong Zhang /* members in merge */ 418445158ffSHong Zhang ierr = PetscFree(id_r);CHKERRQ(ierr); 419445158ffSHong Zhang ierr = PetscFree(len_s);CHKERRQ(ierr); 420445158ffSHong Zhang ierr = PetscFree(len_r);CHKERRQ(ierr); 421445158ffSHong Zhang ierr = PetscFree(buf_ri[0]);CHKERRQ(ierr); 422445158ffSHong Zhang ierr = PetscFree(buf_ri);CHKERRQ(ierr); 423445158ffSHong Zhang ierr = PetscFree(buf_rj[0]);CHKERRQ(ierr); 424445158ffSHong Zhang ierr = PetscFree(buf_rj);CHKERRQ(ierr); 425445158ffSHong Zhang ierr = PetscLayoutDestroy(&rowmap);CHKERRQ(ierr); 42615a3b8e2SHong Zhang 42715a3b8e2SHong Zhang /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */ 42815a3b8e2SHong Zhang Cmpi->assembled = PETSC_FALSE; 42915a3b8e2SHong Zhang Cmpi->ops->destroy = MatDestroy_MPIAIJ_PtAP; 430445158ffSHong Zhang Cmpi->ops->duplicate = 0; //MatDuplicate_MPIAIJ_MatPtAP; //not done yet 431748c7196SHong Zhang Cmpi->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ_new; 43215a3b8e2SHong Zhang 43315a3b8e2SHong Zhang /* attach the supporting struct to Cmpi for reuse */ 43415a3b8e2SHong Zhang c = (Mat_MPIAIJ*)Cmpi->data; 43515a3b8e2SHong Zhang c->ptap = ptap; 436445158ffSHong Zhang ptap->api = api; 437445158ffSHong Zhang ptap->apj = apj; 438445158ffSHong Zhang ptap->apv = apv; 43915a3b8e2SHong Zhang ptap->rmax = ap_rmax; 44015a3b8e2SHong Zhang *C = Cmpi; 44115a3b8e2SHong Zhang 4422259aa2eSHong Zhang /* flag 'scalable' determines which implementations to be used: 4432259aa2eSHong Zhang 0: do dense axpy in MatPtAPNumeric() - fast, but requires storage of a nonscalable dense array apa; 4442259aa2eSHong Zhang 1: do sparse axpy in MatPtAPNumeric() - might slow, uses a sparse array apa */ 4452259aa2eSHong Zhang /* set default scalable */ 4462259aa2eSHong Zhang ptap->scalable = PETSC_FALSE; //PETSC_TRUE; 4472259aa2eSHong Zhang 44815a3b8e2SHong Zhang ierr = PetscOptionsGetBool(((PetscObject)(*C))->prefix,"-matptap_scalable",&ptap->scalable,NULL);CHKERRQ(ierr); 4492259aa2eSHong Zhang if (!ptap->scalable) { /* Do dense axpy */ 4502259aa2eSHong Zhang ierr = PetscCalloc1(P->cmap->N,&ptap->apa);CHKERRQ(ierr); 4512259aa2eSHong Zhang } else { 4522259aa2eSHong Zhang //ierr = PetscCalloc1(ap_rmax+1,&ptap->apa);CHKERRQ(ierr); 4532259aa2eSHong Zhang } 4542259aa2eSHong Zhang 45580bb4639SHong Zhang ierr = PetscTime(&t4);CHKERRQ(ierr); 456e9c1f85fSHong Zhang if (rank == 1) { 457445158ffSHong Zhang printf("PtAPSym: %g + %g + %g + %g + %g + %g = %g\n",t11-t0,t1-t11,t12-t11,t2-t2,t3-t2,t4-t3,t4-t0); 458e9c1f85fSHong Zhang } 4590d3441aeSHong Zhang PetscFunctionReturn(0); 4600d3441aeSHong Zhang } 4610d3441aeSHong Zhang 4620d3441aeSHong Zhang #undef __FUNCT__ 4630d3441aeSHong Zhang #define __FUNCT__ "MatPtAPNumeric_MPIAIJ_MPIAIJ_new" 4640d3441aeSHong Zhang PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ_new(Mat A,Mat P,Mat C) 4650d3441aeSHong Zhang { 4660d3441aeSHong Zhang PetscErrorCode ierr; 4670d3441aeSHong Zhang Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data; 4680d3441aeSHong Zhang Mat_SeqAIJ *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data; 4690d3441aeSHong Zhang Mat_PtAPMPI *ptap = c->ptap; 4709ce11a7cSHong Zhang Mat AP_loc,C_loc,C_oth; 4710d3441aeSHong Zhang PetscInt i,rstart,rend,cm,ncols,row; 4720d3441aeSHong Zhang PetscMPIInt rank; 4730d3441aeSHong Zhang MPI_Comm comm; 4740d3441aeSHong Zhang const PetscInt *cols; 4750d3441aeSHong Zhang const PetscScalar *vals; 4760d3441aeSHong Zhang PetscLogDouble t0,t1,t2,t3,t4,eR,eAP,eCseq,eCmpi; 4770d3441aeSHong Zhang 4780d3441aeSHong Zhang PetscFunctionBegin; 4790d3441aeSHong Zhang ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr); 4800d3441aeSHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 4810d3441aeSHong Zhang 4820d3441aeSHong Zhang ierr = MatZeroEntries(C);CHKERRQ(ierr); 4830d3441aeSHong Zhang 484e497d3c8SHong Zhang /* 1) get R = Pd^T,Ro = Po^T */ 4850d3441aeSHong Zhang ierr = PetscTime(&t0);CHKERRQ(ierr); 486748c7196SHong Zhang if (ptap->reuse == MAT_REUSE_MATRIX) { 4870d3441aeSHong Zhang ierr = MatTranspose_SeqAIJ(p->A,MAT_REUSE_MATRIX,&ptap->Rd);CHKERRQ(ierr); 4880d3441aeSHong Zhang ierr = MatTranspose_SeqAIJ(p->B,MAT_REUSE_MATRIX,&ptap->Ro);CHKERRQ(ierr); 489748c7196SHong Zhang } 4900d3441aeSHong Zhang ierr = PetscTime(&t1);CHKERRQ(ierr); 4910d3441aeSHong Zhang eR = t1 - t0; 4920d3441aeSHong Zhang 493e497d3c8SHong Zhang /* 2) get AP_loc */ 4940d3441aeSHong Zhang AP_loc = ptap->AP_loc; 495e497d3c8SHong Zhang Mat_SeqAIJ *ap=(Mat_SeqAIJ*)AP_loc->data; 4960d3441aeSHong Zhang 497e497d3c8SHong Zhang /* 2-1) get P_oth = ptap->P_oth and P_loc = ptap->P_loc */ 4980d3441aeSHong Zhang /*-----------------------------------------------------*/ 499748c7196SHong Zhang if (ptap->reuse == MAT_REUSE_MATRIX) { 500748c7196SHong Zhang /* P_oth and P_loc are obtained in MatPtASymbolic() when reuse == MAT_INITIAL_MATRIX */ 5010d3441aeSHong Zhang ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 5020d3441aeSHong Zhang ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 503e497d3c8SHong Zhang } 5040d3441aeSHong Zhang 505e497d3c8SHong Zhang /* 2-2) compute numeric C_seq = P_loc^T*A_loc*P - dominating part */ 5060d3441aeSHong Zhang /*--------------------------------------------------------------*/ 5070d3441aeSHong Zhang /* get data from symbolic products */ 5080d3441aeSHong Zhang Mat_SeqAIJ *p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 5090d3441aeSHong Zhang Mat_SeqAIJ *p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 510e497d3c8SHong Zhang PetscInt *api,*apj,am = A->rmap->n,j,col,apnz; 511e497d3c8SHong Zhang PetscScalar *apa = ptap->apa; 512445158ffSHong Zhang //if (ptap->reuse == MAT_REUSE_MATRIX) { 513445158ffSHong Zhang api = ptap->api; //ap->i; 514445158ffSHong Zhang apj = ptap->apj; //ap->j; 515e497d3c8SHong Zhang for (i=0; i<am; i++) { 516445158ffSHong Zhang /* AP[i,:] = A[i,:]*P = Ad*P_loc Ao*P_oth */ 517e497d3c8SHong Zhang AProw_nonscalable(i,ad,ao,p_loc,p_oth,apa); 518e497d3c8SHong Zhang apnz = api[i+1] - api[i]; 519e497d3c8SHong Zhang for (j=0; j<apnz; j++) { 520e497d3c8SHong Zhang col = apj[j+api[i]]; 521e497d3c8SHong Zhang ap->a[j+ap->i[i]] = apa[col]; 522e497d3c8SHong Zhang apa[col] = 0.0; 5230d3441aeSHong Zhang } 5240d3441aeSHong Zhang } 525445158ffSHong Zhang //} 5260d3441aeSHong Zhang 5270d3441aeSHong Zhang ierr = PetscTime(&t2);CHKERRQ(ierr); 5280d3441aeSHong Zhang eAP = t2 - t1; 5290d3441aeSHong Zhang 530e497d3c8SHong Zhang /* 3) C_loc = R*AP_loc, Co = Ro*AP_loc */ 531445158ffSHong Zhang ierr = ((ptap->C_loc)->ops->matmultnumeric)(ptap->Rd,AP_loc,ptap->C_loc);CHKERRQ(ierr); 532445158ffSHong Zhang ierr = ((ptap->C_oth)->ops->matmultnumeric)(ptap->Ro,AP_loc,ptap->C_oth);CHKERRQ(ierr); 533445158ffSHong Zhang 5349ce11a7cSHong Zhang C_loc = ptap->C_loc; 5359ce11a7cSHong Zhang C_oth = ptap->C_oth; 5360d3441aeSHong Zhang //printf("[%d] Co %d, %d\n", rank,Co->rmap->N,Co->cmap->N); 5370d3441aeSHong Zhang ierr = PetscTime(&t3);CHKERRQ(ierr); 5380d3441aeSHong Zhang eCseq = t3 - t2; 5390d3441aeSHong Zhang 5400d3441aeSHong Zhang /* add C_loc and Co to to C */ 5410d3441aeSHong Zhang ierr = MatGetOwnershipRange(C,&rstart,&rend);CHKERRQ(ierr); 5420d3441aeSHong Zhang 5430d3441aeSHong Zhang /* C_loc -> C */ 5440d3441aeSHong Zhang cm = C_loc->rmap->N; 5459ce11a7cSHong Zhang Mat_SeqAIJ *c_seq; 5469ce11a7cSHong Zhang c_seq = (Mat_SeqAIJ*)C_loc->data; 5470d3441aeSHong Zhang for (i=0; i<cm; i++) { 5489ce11a7cSHong Zhang ncols = c_seq->i[i+1] - c_seq->i[i]; 5490d3441aeSHong Zhang row = rstart + i; 5509ce11a7cSHong Zhang cols = c_seq->j + c_seq->i[i]; 5519ce11a7cSHong Zhang vals = c_seq->a + c_seq->i[i]; 5520d3441aeSHong Zhang ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr); 5530d3441aeSHong Zhang } 5540d3441aeSHong Zhang 5550d3441aeSHong Zhang /* Co -> C, off-processor part */ 5560d3441aeSHong Zhang //printf("[%d] p->B %d, %d\n",rank,p->B->rmap->N,p->B->cmap->N); 5579ce11a7cSHong Zhang cm = C_oth->rmap->N; 5589ce11a7cSHong Zhang c_seq = (Mat_SeqAIJ*)C_oth->data; 5599ce11a7cSHong Zhang for (i=0; i<cm; i++) { 5609ce11a7cSHong Zhang ncols = c_seq->i[i+1] - c_seq->i[i]; 5610d3441aeSHong Zhang row = p->garray[i]; 5629ce11a7cSHong Zhang cols = c_seq->j + c_seq->i[i]; 5639ce11a7cSHong Zhang vals = c_seq->a + c_seq->i[i]; 5640d3441aeSHong Zhang //printf("[%d] row[%d] = %d\n",rank,i,row); 5650d3441aeSHong Zhang ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr); 5660d3441aeSHong Zhang } 5670d3441aeSHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 5680d3441aeSHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 5690d3441aeSHong Zhang ierr = PetscTime(&t4);CHKERRQ(ierr); 5700d3441aeSHong Zhang eCmpi = t4 - t3; 5710d3441aeSHong Zhang 5720d3441aeSHong Zhang if (rank==1) { 5730d3441aeSHong Zhang ierr = PetscPrintf(MPI_COMM_SELF," R %g, AP %g, Cseq %g, Cmpi %g = %g\n", eR,eAP,eCseq,eCmpi,eR+eAP+eCseq+eCmpi);CHKERRQ(ierr); 5740d3441aeSHong Zhang } 575748c7196SHong Zhang ptap->reuse = MAT_REUSE_MATRIX; 5760d3441aeSHong Zhang PetscFunctionReturn(0); 5770d3441aeSHong Zhang } 5780d3441aeSHong Zhang 5790d3441aeSHong Zhang #undef __FUNCT__ 58042c57489SHong Zhang #define __FUNCT__ "MatPtAPSymbolic_MPIAIJ_MPIAIJ" 58142c57489SHong Zhang PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C) 58242c57489SHong Zhang { 58342c57489SHong Zhang PetscErrorCode ierr; 58477584fe6SHong Zhang Mat Cmpi; 585f8487c73SHong Zhang Mat_PtAPMPI *ptap; 5860298fd71SBarry Smith PetscFreeSpaceList free_space=NULL,current_space=NULL; 587f8487c73SHong Zhang Mat_MPIAIJ *a =(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c; 58842c57489SHong Zhang Mat_SeqAIJ *ad =(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data; 58942c57489SHong Zhang Mat_SeqAIJ *p_loc,*p_oth; 590ded4f561SHong Zhang PetscInt *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pdti,*pdtj,*poti,*potj,*ptJ; 59177584fe6SHong Zhang PetscInt *adi=ad->i,*aj,*aoi=ao->i,nnz; 592a914927fSHong Zhang PetscInt *lnk,*owners_co,*coi,*coj,i,k,pnz,row; 593d16ca5f0SHong Zhang PetscInt am=A->rmap->n,pN=P->cmap->N,pm=P->rmap->n,pn=P->cmap->n; 59442c57489SHong Zhang PetscBT lnkbt; 595ce94432eSBarry Smith MPI_Comm comm; 596a914927fSHong Zhang PetscMPIInt size,rank,tagi,tagj,*len_si,*len_s,*len_ri,icompleted=0; 59742c57489SHong Zhang PetscInt **buf_rj,**buf_ri,**buf_ri_k; 59842c57489SHong Zhang PetscInt len,proc,*dnz,*onz,*owners; 59924ecddacSHong Zhang PetscInt nzi,*pti,*ptj; 60042c57489SHong Zhang PetscInt nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci; 601ded4f561SHong Zhang MPI_Request *swaits,*rwaits; 602ded4f561SHong Zhang MPI_Status *sstatus,rstatus; 60342c57489SHong Zhang Mat_Merge_SeqsToMPI *merge; 60477584fe6SHong Zhang PetscInt *api,*apj,*Jptr,apnz,*prmap=p->garray,pon,nspacedouble=0,j,ap_rmax=0; 605d16ca5f0SHong Zhang PetscReal afill=1.0,afill_tmp; 606a914927fSHong Zhang PetscInt rmax; 60742c57489SHong Zhang 60842c57489SHong Zhang PetscFunctionBegin; 609ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 61024ecddacSHong Zhang 611c5af039cSHong Zhang /* check if matrix local sizes are compatible */ 612c5af039cSHong Zhang if (A->rmap->rstart != P->rmap->rstart || A->rmap->rend != P->rmap->rend) { 613c5af039cSHong Zhang SETERRQ4(comm,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, Arow (%D, %D) != Prow (%D,%D)",A->rmap->rstart,A->rmap->rend,P->rmap->rstart,P->rmap->rend); 614c5af039cSHong Zhang } 615c5af039cSHong Zhang if (A->cmap->rstart != P->rmap->rstart || A->cmap->rend != P->rmap->rend) { 616c5af039cSHong Zhang SETERRQ4(comm,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, Acol (%D, %D) != Prow (%D,%D)",A->cmap->rstart,A->cmap->rend,P->rmap->rstart,P->rmap->rend); 617c5af039cSHong Zhang } 618c5af039cSHong Zhang 61983445d95SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 62083445d95SHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 62183445d95SHong Zhang 622f8487c73SHong Zhang /* create struct Mat_PtAPMPI and attached it to C later */ 623b00a9115SJed Brown ierr = PetscNew(&ptap);CHKERRQ(ierr); 624b00a9115SJed Brown ierr = PetscNew(&merge);CHKERRQ(ierr); 6259f4155fbSHong Zhang ptap->merge = merge; 626f8487c73SHong Zhang ptap->reuse = MAT_INITIAL_MATRIX; 6276c6e00e0SHong Zhang 6286c6e00e0SHong Zhang /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */ 629b7f45c76SHong Zhang ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 630fe615159SHong Zhang 6316c6e00e0SHong Zhang /* get P_loc by taking all local rows of P */ 632a1a4e84aSHong Zhang ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 6336c6e00e0SHong Zhang 634a1a4e84aSHong Zhang p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 635a1a4e84aSHong Zhang p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 6366c6e00e0SHong Zhang pi_loc = p_loc->i; pj_loc = p_loc->j; 6376c6e00e0SHong Zhang pi_oth = p_oth->i; pj_oth = p_oth->j; 63842c57489SHong Zhang 6392addb229SHong Zhang /* (1) compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth (api,apj) */ 6402addb229SHong Zhang /*--------------------------------------------------------------------------*/ 641854ce69bSBarry Smith ierr = PetscMalloc1(am+1,&api);CHKERRQ(ierr); 64282412ba7SHong Zhang api[0] = 0; 64383445d95SHong Zhang 64483445d95SHong Zhang /* create and initialize a linked list */ 645a914927fSHong Zhang ierr = PetscLLCondensedCreate(pN,pN,&lnk,&lnkbt);CHKERRQ(ierr); 64683445d95SHong Zhang 647a914927fSHong Zhang /* Initial FreeSpace size is fill*(nnz(A) + nnz(P)) -OOM for ex56, np=8k on Intrepid! */ 648d16ca5f0SHong Zhang ierr = PetscFreeSpaceGet((PetscInt)(fill*(adi[am]+aoi[am]+pi_loc[pm])),&free_space);CHKERRQ(ierr); 649f4a743e1SHong Zhang current_space = free_space; 650f4a743e1SHong Zhang 651f4a743e1SHong Zhang for (i=0; i<am; i++) { 652f4a743e1SHong Zhang /* diagonal portion of A */ 653ded4f561SHong Zhang nzi = adi[i+1] - adi[i]; 65477584fe6SHong Zhang aj = ad->j + adi[i]; 655ded4f561SHong Zhang for (j=0; j<nzi; j++) { 65677584fe6SHong Zhang row = aj[j]; 65782412ba7SHong Zhang pnz = pi_loc[row+1] - pi_loc[row]; 658ded4f561SHong Zhang Jptr = pj_loc + pi_loc[row]; 65983445d95SHong Zhang /* add non-zero cols of P into the sorted linked list lnk */ 660a914927fSHong Zhang ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 661f4a743e1SHong Zhang } 662f4a743e1SHong Zhang /* off-diagonal portion of A */ 663ded4f561SHong Zhang nzi = aoi[i+1] - aoi[i]; 66477584fe6SHong Zhang aj = ao->j + aoi[i]; 665ded4f561SHong Zhang for (j=0; j<nzi; j++) { 66677584fe6SHong Zhang row = aj[j]; 66782412ba7SHong Zhang pnz = pi_oth[row+1] - pi_oth[row]; 668ded4f561SHong Zhang Jptr = pj_oth + pi_oth[row]; 669a914927fSHong Zhang ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 670f4a743e1SHong Zhang } 671a914927fSHong Zhang apnz = lnk[0]; 67282412ba7SHong Zhang api[i+1] = api[i] + apnz; 67377584fe6SHong Zhang if (ap_rmax < apnz) ap_rmax = apnz; 674f4a743e1SHong Zhang 67583445d95SHong Zhang /* if free space is not available, double the total space in the list */ 67682412ba7SHong Zhang if (current_space->local_remaining<apnz) { 6772ba1acd5SHong Zhang ierr = PetscFreeSpaceGet(apnz+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 678f2b054eeSHong Zhang nspacedouble++; 679f4a743e1SHong Zhang } 680f4a743e1SHong Zhang 681f4a743e1SHong Zhang /* Copy data into free space, then initialize lnk */ 682a914927fSHong Zhang ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 6832205254eSKarl Rupp 68482412ba7SHong Zhang current_space->array += apnz; 68582412ba7SHong Zhang current_space->local_used += apnz; 68682412ba7SHong Zhang current_space->local_remaining -= apnz; 687f4a743e1SHong Zhang } 688a914927fSHong Zhang 68982412ba7SHong Zhang /* Allocate space for apj, initialize apj, and */ 690f4a743e1SHong Zhang /* destroy list of free space and other temporary array(s) */ 691854ce69bSBarry Smith ierr = PetscMalloc1(api[am]+1,&apj);CHKERRQ(ierr); 69277584fe6SHong Zhang ierr = PetscFreeSpaceContiguous(&free_space,apj);CHKERRQ(ierr); 693118727c9SMark F. Adams afill_tmp = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1); 694d16ca5f0SHong Zhang if (afill_tmp > afill) afill = afill_tmp; 695f4a743e1SHong Zhang 6962addb229SHong Zhang /* (2) determine symbolic Co=(p->B)^T*AP - send to others (coi,coj)*/ 6972addb229SHong Zhang /*-----------------------------------------------------------------*/ 698de0260b3SHong Zhang ierr = MatGetSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr); 69942c57489SHong Zhang 700ded4f561SHong Zhang /* then, compute symbolic Co = (p->B)^T*AP */ 701d0f46423SBarry Smith pon = (p->B)->cmap->n; /* total num of rows to be sent to other processors 70283445d95SHong Zhang >= (num of nonzero rows of C_seq) - pn */ 703854ce69bSBarry Smith ierr = PetscMalloc1(pon+1,&coi);CHKERRQ(ierr); 704de0260b3SHong Zhang coi[0] = 0; 70542c57489SHong Zhang 706d16ca5f0SHong Zhang /* set initial free space to be fill*(nnz(p->B) + nnz(AP)) */ 707d16ca5f0SHong Zhang nnz = fill*(poti[pon] + api[am]); 70822d28d08SBarry Smith ierr = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr); 70942c57489SHong Zhang current_space = free_space; 71042c57489SHong Zhang 711de0260b3SHong Zhang for (i=0; i<pon; i++) { 71282412ba7SHong Zhang pnz = poti[i+1] - poti[i]; 71377584fe6SHong Zhang ptJ = potj + poti[i]; 71477584fe6SHong Zhang for (j=0; j<pnz; j++) { 71577584fe6SHong Zhang row = ptJ[j]; /* row of AP == col of Pot */ 71682412ba7SHong Zhang apnz = api[row+1] - api[row]; 717ded4f561SHong Zhang Jptr = apj + api[row]; 71883445d95SHong Zhang /* add non-zero cols of AP into the sorted linked list lnk */ 719a914927fSHong Zhang ierr = PetscLLCondensedAddSorted(apnz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 72042c57489SHong Zhang } 721a914927fSHong Zhang nnz = lnk[0]; 72242c57489SHong Zhang 72383445d95SHong Zhang /* If free space is not available, double the total space in the list */ 724ded4f561SHong Zhang if (current_space->local_remaining<nnz) { 7254238b7adSHong Zhang ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 726d16ca5f0SHong Zhang nspacedouble++; 72742c57489SHong Zhang } 72842c57489SHong Zhang 72942c57489SHong Zhang /* Copy data into free space, and zero out denserows */ 730a914927fSHong Zhang ierr = PetscLLCondensedClean(pN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 7312205254eSKarl Rupp 732ded4f561SHong Zhang current_space->array += nnz; 733ded4f561SHong Zhang current_space->local_used += nnz; 734ded4f561SHong Zhang current_space->local_remaining -= nnz; 7352205254eSKarl Rupp 736ded4f561SHong Zhang coi[i+1] = coi[i] + nnz; 73742c57489SHong Zhang } 7386b911d16SHong Zhang 7396b911d16SHong Zhang ierr = PetscMalloc1(coi[pon],&coj);CHKERRQ(ierr); 740a1a86e44SBarry Smith ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr); 741118727c9SMark F. Adams afill_tmp = (PetscReal)coi[pon]/(poti[pon] + api[am]+1); 742d16ca5f0SHong Zhang if (afill_tmp > afill) afill = afill_tmp; 743de0260b3SHong Zhang ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr); 74442c57489SHong Zhang 7456b911d16SHong Zhang /* (3) send j-array (coj) of Co to other processors */ 7466b911d16SHong Zhang /*--------------------------------------------------*/ 7476b911d16SHong Zhang ierr = PetscCalloc1(size,&merge->len_s);CHKERRQ(ierr); 7486b911d16SHong Zhang len_s = merge->len_s; 7496b911d16SHong Zhang merge->nsend = 0; 7506b911d16SHong Zhang 7516b911d16SHong Zhang 75242c57489SHong Zhang /* determine row ownership */ 75326283091SBarry Smith ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr); 7547a2fc3feSBarry Smith merge->rowmap->n = pn; 7557a2fc3feSBarry Smith merge->rowmap->bs = 1; 7562205254eSKarl Rupp 75726283091SBarry Smith ierr = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr); 7587a2fc3feSBarry Smith owners = merge->rowmap->range; 75942c57489SHong Zhang 76042c57489SHong Zhang /* determine the number of messages to send, their lengths */ 761dcca6d9dSJed Brown ierr = PetscMalloc2(size,&len_si,size,&sstatus);CHKERRQ(ierr); 76283445d95SHong Zhang ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 763854ce69bSBarry Smith ierr = PetscMalloc1(size+2,&owners_co);CHKERRQ(ierr); 764de0260b3SHong Zhang 76583445d95SHong Zhang proc = 0; 766de0260b3SHong Zhang for (i=0; i<pon; i++) { 767de0260b3SHong Zhang while (prmap[i] >= owners[proc+1]) proc++; 7682addb229SHong Zhang len_si[proc]++; /* num of rows in Co(=Pt*AP) to be sent to [proc] */ 769ee6e4b5bSHong Zhang len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */ 77083445d95SHong Zhang } 771de0260b3SHong Zhang 7722addb229SHong Zhang len = 0; /* max length of buf_si[], see (4) */ 773de0260b3SHong Zhang owners_co[0] = 0; 77442c57489SHong Zhang for (proc=0; proc<size; proc++) { 775de0260b3SHong Zhang owners_co[proc+1] = owners_co[proc] + len_si[proc]; 776ee6e4b5bSHong Zhang if (len_s[proc]) { 77742c57489SHong Zhang merge->nsend++; 7782addb229SHong Zhang len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */ 77942c57489SHong Zhang len += len_si[proc]; 78042c57489SHong Zhang } 78142c57489SHong Zhang } 78242c57489SHong Zhang 783ded4f561SHong Zhang /* determine the number and length of messages to receive for coi and coj */ 7840298fd71SBarry Smith ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&merge->nrecv);CHKERRQ(ierr); 78542c57489SHong Zhang ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr); 78642c57489SHong Zhang 787ded4f561SHong Zhang /* post the Irecv and Isend of coj */ 788529bc5dcSHong Zhang ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr); 789529bc5dcSHong Zhang ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr); 790854ce69bSBarry Smith ierr = PetscMalloc1(merge->nsend+1,&swaits);CHKERRQ(ierr); 79142c57489SHong Zhang for (proc=0, k=0; proc<size; proc++) { 79242c57489SHong Zhang if (!len_s[proc]) continue; 793de0260b3SHong Zhang i = owners_co[proc]; 794529bc5dcSHong Zhang ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr); 79542c57489SHong Zhang k++; 79642c57489SHong Zhang } 79742c57489SHong Zhang 798ded4f561SHong Zhang /* receives and sends of coj are complete */ 79977584fe6SHong Zhang for (i=0; i<merge->nrecv; i++) { 800c280ed6aSJed Brown ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 801ded4f561SHong Zhang } 802ded4f561SHong Zhang ierr = PetscFree(rwaits);CHKERRQ(ierr); 8030c468ba9SBarry Smith if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);} 80482412ba7SHong Zhang 8056b911d16SHong Zhang /* (4) send and recv coi */ 8066b911d16SHong Zhang /*-----------------------*/ 807529bc5dcSHong Zhang ierr = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr); 808529bc5dcSHong Zhang ierr = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr); 809854ce69bSBarry Smith ierr = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr); 81042c57489SHong Zhang buf_si = buf_s; /* points to the beginning of k-th msg to be sent */ 81142c57489SHong Zhang for (proc=0,k=0; proc<size; proc++) { 81242c57489SHong Zhang if (!len_s[proc]) continue; 81342c57489SHong Zhang /* form outgoing message for i-structure: 81442c57489SHong Zhang buf_si[0]: nrows to be sent 81542c57489SHong Zhang [1:nrows]: row index (global) 81642c57489SHong Zhang [nrows+1:2*nrows+1]: i-structure index 81742c57489SHong Zhang */ 81842c57489SHong Zhang /*-------------------------------------------*/ 8192addb229SHong Zhang nrows = len_si[proc]/2 - 1; /* num of rows in Co to be sent to [proc] */ 82042c57489SHong Zhang buf_si_i = buf_si + nrows+1; 82142c57489SHong Zhang buf_si[0] = nrows; 82242c57489SHong Zhang buf_si_i[0] = 0; 82342c57489SHong Zhang nrows = 0; 824de0260b3SHong Zhang for (i=owners_co[proc]; i<owners_co[proc+1]; i++) { 825ded4f561SHong Zhang nzi = coi[i+1] - coi[i]; 826ded4f561SHong Zhang buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */ 827de0260b3SHong Zhang buf_si[nrows+1] = prmap[i] -owners[proc]; /* local row index */ 82842c57489SHong Zhang nrows++; 82942c57489SHong Zhang } 830529bc5dcSHong Zhang ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr); 83142c57489SHong Zhang k++; 83242c57489SHong Zhang buf_si += len_si[proc]; 83342c57489SHong Zhang } 834ded4f561SHong Zhang i = merge->nrecv; 835ded4f561SHong Zhang while (i--) { 836c280ed6aSJed Brown ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 837ded4f561SHong Zhang } 838ded4f561SHong Zhang ierr = PetscFree(rwaits);CHKERRQ(ierr); 8390c468ba9SBarry Smith if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);} 840a914927fSHong Zhang 84124ecddacSHong Zhang ierr = PetscFree2(len_si,sstatus);CHKERRQ(ierr); 84242c57489SHong Zhang ierr = PetscFree(len_ri);CHKERRQ(ierr); 843ded4f561SHong Zhang ierr = PetscFree(swaits);CHKERRQ(ierr); 84442c57489SHong Zhang ierr = PetscFree(buf_s);CHKERRQ(ierr); 84542c57489SHong Zhang 8466b911d16SHong Zhang /* (5) compute the local portion of C (mpi mat) */ 8476b911d16SHong Zhang /*----------------------------------------------*/ 848ded4f561SHong Zhang ierr = MatGetSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr); 849ded4f561SHong Zhang 85024ecddacSHong Zhang /* allocate pti array and free space for accumulating nonzero column info */ 851854ce69bSBarry Smith ierr = PetscMalloc1(pn+1,&pti);CHKERRQ(ierr); 85224ecddacSHong Zhang pti[0] = 0; 85342c57489SHong Zhang 854d16ca5f0SHong Zhang /* set initial free space to be fill*(nnz(P) + nnz(AP)) */ 855d16ca5f0SHong Zhang nnz = fill*(pi_loc[pm] + api[am]); 85622d28d08SBarry Smith ierr = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr); 85742c57489SHong Zhang current_space = free_space; 85842c57489SHong Zhang 859dcca6d9dSJed Brown ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr); 86042c57489SHong Zhang for (k=0; k<merge->nrecv; k++) { 86142c57489SHong Zhang buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 86242c57489SHong Zhang nrows = *buf_ri_k[k]; 86342c57489SHong Zhang nextrow[k] = buf_ri_k[k] + 1; /* next row number of k-th recved i-structure */ 86442c57489SHong Zhang nextci[k] = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure */ 86542c57489SHong Zhang } 86642c57489SHong Zhang ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr); 86708cb4532SHong Zhang rmax = 0; 86842c57489SHong Zhang for (i=0; i<pn; i++) { 869ded4f561SHong Zhang /* add pdt[i,:]*AP into lnk */ 870ded4f561SHong Zhang pnz = pdti[i+1] - pdti[i]; 87177584fe6SHong Zhang ptJ = pdtj + pdti[i]; 87277584fe6SHong Zhang for (j=0; j<pnz; j++) { 87377584fe6SHong Zhang row = ptJ[j]; /* row of AP == col of Pt */ 874ded4f561SHong Zhang apnz = api[row+1] - api[row]; 875ded4f561SHong Zhang Jptr = apj + api[row]; 876ded4f561SHong Zhang /* add non-zero cols of AP into the sorted linked list lnk */ 877a914927fSHong Zhang ierr = PetscLLCondensedAddSorted(apnz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 878ded4f561SHong Zhang } 879a914927fSHong Zhang 88042c57489SHong Zhang /* add received col data into lnk */ 88142c57489SHong Zhang for (k=0; k<merge->nrecv; k++) { /* k-th received message */ 88242c57489SHong Zhang if (i == *nextrow[k]) { /* i-th row */ 883ded4f561SHong Zhang nzi = *(nextci[k]+1) - *nextci[k]; 884ded4f561SHong Zhang Jptr = buf_rj[k] + *nextci[k]; 885a914927fSHong Zhang ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr); 88642c57489SHong Zhang nextrow[k]++; nextci[k]++; 88742c57489SHong Zhang } 88842c57489SHong Zhang } 889a914927fSHong Zhang nnz = lnk[0]; 89042c57489SHong Zhang 89142c57489SHong Zhang /* if free space is not available, make more free space */ 892ded4f561SHong Zhang if (current_space->local_remaining<nnz) { 8934238b7adSHong Zhang ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 894d16ca5f0SHong Zhang nspacedouble++; 89542c57489SHong Zhang } 89642c57489SHong Zhang /* copy data into free space, then initialize lnk */ 897a914927fSHong Zhang ierr = PetscLLCondensedClean(pN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 898ded4f561SHong Zhang ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr); 8992205254eSKarl Rupp 900ded4f561SHong Zhang current_space->array += nnz; 901ded4f561SHong Zhang current_space->local_used += nnz; 902ded4f561SHong Zhang current_space->local_remaining -= nnz; 9032205254eSKarl Rupp 90424ecddacSHong Zhang pti[i+1] = pti[i] + nnz; 90508cb4532SHong Zhang if (nnz > rmax) rmax = nnz; 90642c57489SHong Zhang } 907ded4f561SHong Zhang ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr); 9080572522cSBarry Smith ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr); 90942c57489SHong Zhang 910854ce69bSBarry Smith ierr = PetscMalloc1(pti[pn]+1,&ptj);CHKERRQ(ierr); 91124ecddacSHong Zhang ierr = PetscFreeSpaceContiguous(&free_space,ptj);CHKERRQ(ierr); 91224ecddacSHong Zhang afill_tmp = (PetscReal)pti[pn]/(pi_loc[pm] + api[am]+1); 913d16ca5f0SHong Zhang if (afill_tmp > afill) afill = afill_tmp; 91442c57489SHong Zhang ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 91542c57489SHong Zhang 9166b911d16SHong Zhang /* (6) create symbolic parallel matrix Cmpi */ 9176b911d16SHong Zhang /*------------------------------------------*/ 91877584fe6SHong Zhang ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr); 91977584fe6SHong Zhang ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 92033d57670SJed Brown ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr); 92177584fe6SHong Zhang ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr); 92277584fe6SHong Zhang ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr); 92342c57489SHong Zhang ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 924a2f3521dSMark F. Adams 925b091e509SHong Zhang merge->bi = pti; /* Cseq->i */ 926b091e509SHong Zhang merge->bj = ptj; /* Cseq->j */ 927b091e509SHong Zhang merge->coi = coi; /* Co->i */ 928b091e509SHong Zhang merge->coj = coj; /* Co->j */ 92942c57489SHong Zhang merge->buf_ri = buf_ri; 93042c57489SHong Zhang merge->buf_rj = buf_rj; 931de0260b3SHong Zhang merge->owners_co = owners_co; 93277584fe6SHong Zhang merge->destroy = Cmpi->ops->destroy; 93377584fe6SHong Zhang merge->duplicate = Cmpi->ops->duplicate; 934cc6cb787SHong Zhang 93577584fe6SHong Zhang /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */ 936a914927fSHong Zhang Cmpi->assembled = PETSC_FALSE; 93777584fe6SHong Zhang Cmpi->ops->destroy = MatDestroy_MPIAIJ_PtAP; 93877584fe6SHong Zhang Cmpi->ops->duplicate = MatDuplicate_MPIAIJ_MatPtAP; 93942c57489SHong Zhang 94077584fe6SHong Zhang /* attach the supporting struct to Cmpi for reuse */ 94177584fe6SHong Zhang c = (Mat_MPIAIJ*)Cmpi->data; 942f8487c73SHong Zhang c->ptap = ptap; 94377584fe6SHong Zhang ptap->api = api; 94477584fe6SHong Zhang ptap->apj = apj; 945d6ab1dc0SHong Zhang ptap->rmax = ap_rmax; 94677584fe6SHong Zhang *C = Cmpi; 947414894bdSHong Zhang 948414894bdSHong Zhang /* flag 'scalable' determines which implementations to be used: 949414894bdSHong Zhang 0: do dense axpy in MatPtAPNumeric() - fast, but requires storage of a nonscalable dense array apa; 950414894bdSHong Zhang 1: do sparse axpy in MatPtAPNumeric() - might slow, uses a sparse array apa */ 951414894bdSHong Zhang /* set default scalable */ 9529516a85cSHong Zhang ptap->scalable = PETSC_FALSE; //PETSC_TRUE; 9532205254eSKarl Rupp 9540298fd71SBarry Smith ierr = PetscOptionsGetBool(((PetscObject)Cmpi)->prefix,"-matptap_scalable",&ptap->scalable,NULL);CHKERRQ(ierr); 955414894bdSHong Zhang if (!ptap->scalable) { /* Do dense axpy */ 9561795a4d1SJed Brown ierr = PetscCalloc1(pN,&ptap->apa);CHKERRQ(ierr); 957414894bdSHong Zhang } else { 9581795a4d1SJed Brown ierr = PetscCalloc1(ap_rmax+1,&ptap->apa);CHKERRQ(ierr); 959414894bdSHong Zhang } 960414894bdSHong Zhang 961f2b054eeSHong Zhang #if defined(PETSC_USE_INFO) 96224ecddacSHong Zhang if (pti[pn] != 0) { 96357622a8eSBarry Smith ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr); 96457622a8eSBarry Smith ierr = PetscInfo1(Cmpi,"Use MatPtAP(A,P,MatReuse,%g,&C) for best performance.\n",(double)afill);CHKERRQ(ierr); 965f2b054eeSHong Zhang } else { 96677584fe6SHong Zhang ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr); 967f2b054eeSHong Zhang } 968f2b054eeSHong Zhang #endif 96942c57489SHong Zhang PetscFunctionReturn(0); 97042c57489SHong Zhang } 97142c57489SHong Zhang 97242c57489SHong Zhang #undef __FUNCT__ 97342c57489SHong Zhang #define __FUNCT__ "MatPtAPNumeric_MPIAIJ_MPIAIJ" 97442c57489SHong Zhang PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C) 97542c57489SHong Zhang { 97642c57489SHong Zhang PetscErrorCode ierr; 977f8487c73SHong Zhang Mat_MPIAIJ *a =(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data; 97842c57489SHong Zhang Mat_SeqAIJ *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data; 979de0260b3SHong Zhang Mat_SeqAIJ *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data; 98042c57489SHong Zhang Mat_SeqAIJ *p_loc,*p_oth; 981f8487c73SHong Zhang Mat_PtAPMPI *ptap; 9829f4155fbSHong Zhang Mat_Merge_SeqsToMPI *merge; 9831d633602SHong Zhang PetscInt *adi=ad->i,*aoi=ao->i,*adj,*aoj,*apJ,nextp; 98482412ba7SHong Zhang PetscInt *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pJ,*pj; 98582412ba7SHong Zhang PetscInt i,j,k,anz,pnz,apnz,nextap,row,*cj; 986e900a757SHong Zhang MatScalar *ada,*aoa,*apa,*pa,*ca,*pa_loc,*pa_oth,valtmp; 987d0f46423SBarry Smith PetscInt am =A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n; 988ce94432eSBarry Smith MPI_Comm comm; 98942c57489SHong Zhang PetscMPIInt size,rank,taga,*len_s; 99082412ba7SHong Zhang PetscInt *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci; 99142c57489SHong Zhang PetscInt **buf_ri,**buf_rj; 99250f5bf86SHong Zhang PetscInt cnz=0,*bj_i,*bi,*bj,bnz,nextcj; /* bi,bj,ba: local array of C(mpi mat) */ 99342c57489SHong Zhang MPI_Request *s_waits,*r_waits; 99442c57489SHong Zhang MPI_Status *status; 99582412ba7SHong Zhang MatScalar **abuf_r,*ba_i,*pA,*coa,*ba; 99682412ba7SHong Zhang PetscInt *api,*apj,*coi,*coj; 997d0f46423SBarry Smith PetscInt *poJ=po->j,*pdJ=pd->j,pcstart=P->cmap->rstart,pcend=P->cmap->rend; 9986ce94e8fSHong Zhang PetscBool scalable; 99938571addSHong Zhang #if defined(PTAP_PROFILE) 1000e497d3c8SHong Zhang PetscLogDouble t0,t1,t2,eP,t3,t4,et2_AP=0.0,ePtAP=0.0,t2_0,t2_1,t2_2; 100138571addSHong Zhang #endif 100242c57489SHong Zhang 100342c57489SHong Zhang PetscFunctionBegin; 1004ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr); 100538571addSHong Zhang #if defined(PTAP_PROFILE) 10068563dfccSBarry Smith ierr = PetscTime(&t0);CHKERRQ(ierr); 100738571addSHong Zhang #endif 100842c57489SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 100942c57489SHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 101042c57489SHong Zhang 1011f8487c73SHong Zhang ptap = c->ptap; 1012ce94432eSBarry Smith if (!ptap) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_ARG_INCOMP,"MatPtAP() has not been called to create matrix C yet, cannot use MAT_REUSE_MATRIX"); 1013f8487c73SHong Zhang merge = ptap->merge; 1014414894bdSHong Zhang apa = ptap->apa; 10156ce94e8fSHong Zhang scalable = ptap->scalable; 10166c6e00e0SHong Zhang 1017a1a4e84aSHong Zhang /* 1) get P_oth = ptap->P_oth and P_loc = ptap->P_loc */ 10186b911d16SHong Zhang /*-----------------------------------------------------*/ 1019f8487c73SHong Zhang if (ptap->reuse == MAT_INITIAL_MATRIX) { 10209f4155fbSHong Zhang /* P_oth and P_loc are obtained in MatPtASymbolic(), skip calling MatGetBrowsOfAoCols() and MatMPIAIJGetLocalMat() */ 1021f8487c73SHong Zhang ptap->reuse = MAT_REUSE_MATRIX; 10226c6e00e0SHong Zhang } else { /* update numerical values of P_oth and P_loc */ 1023b7f45c76SHong Zhang ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 1024a1a4e84aSHong Zhang ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 10256c6e00e0SHong Zhang } 102638571addSHong Zhang #if defined(PTAP_PROFILE) 10278563dfccSBarry Smith ierr = PetscTime(&t1);CHKERRQ(ierr); 102805d62848SHong Zhang eP = t1-t0; 102938571addSHong Zhang #endif 103005d62848SHong Zhang /* 103105d62848SHong Zhang printf("[%d] Ad: %d, %d; Ao: %d, %d; P_loc: %d, %d; P_oth %d, %d;\n",rank, 103205d62848SHong Zhang a->A->rmap->N,a->A->cmap->N,a->B->rmap->N,a->B->cmap->N, 103305d62848SHong Zhang ptap->P_loc->rmap->N,ptap->P_loc->cmap->N, 103405d62848SHong Zhang ptap->P_oth->rmap->N,ptap->P_oth->cmap->N); 103505d62848SHong Zhang */ 1036f8487c73SHong Zhang 1037f9473cf7SHong Zhang /* 2) compute numeric C_seq = P_loc^T*A_loc*P - dominating part */ 1038f9473cf7SHong Zhang /*--------------------------------------------------------------*/ 103942c57489SHong Zhang /* get data from symbolic products */ 1040a1a4e84aSHong Zhang p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 1041a1a4e84aSHong Zhang p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 1042a9d06231SHong Zhang pi_loc=p_loc->i; pj_loc=p_loc->j; pJ=pj_loc; pa_loc=p_loc->a; 104342c57489SHong Zhang pi_oth=p_oth->i; pj_oth=p_oth->j; pa_oth=p_oth->a; 104442c57489SHong Zhang 1045de0260b3SHong Zhang coi = merge->coi; coj = merge->coj; 10461795a4d1SJed Brown ierr = PetscCalloc1(coi[pon]+1,&coa);CHKERRQ(ierr); 1047de0260b3SHong Zhang 1048de0260b3SHong Zhang bi = merge->bi; bj = merge->bj; 10497a2fc3feSBarry Smith owners = merge->rowmap->range; 10501795a4d1SJed Brown ierr = PetscCalloc1(bi[cm]+1,&ba);CHKERRQ(ierr); /* ba: Cseq->a */ 105142c57489SHong Zhang 1052a1a4e84aSHong Zhang api = ptap->api; apj = ptap->apj; 1053d9824c15SHong Zhang 10549516a85cSHong Zhang if (!scalable) { /* Do dense axpy on apa (length of pN, stores A[i,:]*P) - nonscalable, but faster (could take 1/3 scalable time) */ 1055b38d1a2dSJed Brown ierr = PetscInfo(C,"Using non-scalable dense axpy\n");CHKERRQ(ierr); 105605d62848SHong Zhang #if 0 10570d3441aeSHong Zhang /* ------ 10x slower -------------- */ 10580d3441aeSHong Zhang /*==================================*/ 10599516a85cSHong Zhang Mat R = ptap->R; 10609516a85cSHong Zhang Mat_SeqAIJ *r = (Mat_SeqAIJ*)R->data; 10619516a85cSHong Zhang PetscInt *ri=r->i,*rj=r->j,rnz,arow,l,prow,pcol,pN=P->cmap->N; 10629516a85cSHong Zhang PetscScalar *ra=r->a,tmp,cdense[pN]; 10639516a85cSHong Zhang 10649516a85cSHong Zhang ierr = PetscMemzero(cdense,pN*sizeof(PetscScalar));CHKERRQ(ierr); 10659516a85cSHong Zhang for (i=0; i<cm; i++) { /* each row of C or R */ 10669516a85cSHong Zhang rnz = ri[i+1] - ri[i]; 10679516a85cSHong Zhang 10689516a85cSHong Zhang for (j=0; j<rnz; j++) { /* each nz of R */ 10699516a85cSHong Zhang arow = rj[ri[i] + j]; 10709516a85cSHong Zhang 10719516a85cSHong Zhang /* diagonal portion of A */ 10729516a85cSHong Zhang anz = ad->i[arow+1] - ad->i[arow]; 10739516a85cSHong Zhang for (k=0; k<anz; k++) { /* each nz of Ad */ 10749516a85cSHong Zhang tmp = ra[ri[i] + j]*ad->a[ad->i[arow] + k]; 10759516a85cSHong Zhang prow = ad->j[ad->i[arow] + k]; 10769516a85cSHong Zhang pnz = pi_loc[prow+1] - pi_loc[prow]; 10779516a85cSHong Zhang 10789516a85cSHong Zhang for (l=0; l<pnz; l++) { /* each nz of P_loc */ 10799516a85cSHong Zhang pcol = pj_loc[pi_loc[prow] + l]; 10809516a85cSHong Zhang cdense[pcol] += tmp*pa_loc[pi_loc[prow] + l]; 10819516a85cSHong Zhang } 10829516a85cSHong Zhang } 10839516a85cSHong Zhang 10849516a85cSHong Zhang /* off-diagonal portion of A */ 10859516a85cSHong Zhang anz = ao->i[arow+1] - ao->i[arow]; 10869516a85cSHong Zhang for (k=0; k<anz; k++) { /* each nz of Ao */ 10879516a85cSHong Zhang tmp = ra[ri[i] + j]*ao->a[ao->i[arow] + k]; 10889516a85cSHong Zhang prow = ao->j[ao->i[arow] + k]; 10899516a85cSHong Zhang pnz = pi_oth[prow+1] - pi_oth[prow]; 10909516a85cSHong Zhang 10919516a85cSHong Zhang for (l=0; l<pnz; l++) { /* each nz of P_oth */ 10929516a85cSHong Zhang pcol = pj_oth[pi_oth[prow] + l]; 10939516a85cSHong Zhang cdense[pcol] += tmp*pa_oth[pi_oth[prow] + l]; 10949516a85cSHong Zhang } 10959516a85cSHong Zhang } 10969516a85cSHong Zhang 10979516a85cSHong Zhang } //for (j=0; j<rnz; j++) 10989516a85cSHong Zhang 10999516a85cSHong Zhang /* copy cdense[] into ca; zero cdense[] */ 11009516a85cSHong Zhang cnz = bi[i+1] - bi[i]; 11019516a85cSHong Zhang cj = bj + bi[i]; 11029516a85cSHong Zhang ca = ba + bi[i]; 11039516a85cSHong Zhang for (j=0; j<cnz; j++) { 11049516a85cSHong Zhang ca[j] += cdense[cj[j]]; 11059516a85cSHong Zhang cdense[cj[j]] = 0.0; 11069516a85cSHong Zhang } 11079516a85cSHong Zhang #if 0 11089516a85cSHong Zhang if (rank == 0) { 11099516a85cSHong Zhang printf("[%d] row %d: ",rank,i); 11109516a85cSHong Zhang for (j=0; j<pN; j++) printf(" %g,",cdense[j]); 11119516a85cSHong Zhang printf("\n"); 11129516a85cSHong Zhang } 11139516a85cSHong Zhang for (j=0; j<pN; j++) cdense[j]=0.0; // zero cdnese[] 11149516a85cSHong Zhang #endif 11159516a85cSHong Zhang } //for (i=0; i<cm; i++) { 111605d62848SHong Zhang #endif 11179516a85cSHong Zhang 11189516a85cSHong Zhang //========================================== 11199516a85cSHong Zhang 112005d62848SHong Zhang ierr = PetscTime(&t1);CHKERRQ(ierr); 1121a9d06231SHong Zhang for (i=0; i<am; i++) { 1122b091e509SHong Zhang #if defined(PTAP_PROFILE) 11238563dfccSBarry Smith ierr = PetscTime(&t2_0);CHKERRQ(ierr); 1124b091e509SHong Zhang #endif 1125a9d06231SHong Zhang /* 2-a) form i-th sparse row of A_loc*P = Ad*P_loc + Ao*P_oth */ 1126a9d06231SHong Zhang /*------------------------------------------------------------*/ 1127a9d06231SHong Zhang apJ = apj + api[i]; 1128a9d06231SHong Zhang 1129a9d06231SHong Zhang /* diagonal portion of A */ 1130a9d06231SHong Zhang anz = adi[i+1] - adi[i]; 1131a9d06231SHong Zhang adj = ad->j + adi[i]; 1132a9d06231SHong Zhang ada = ad->a + adi[i]; 1133a9d06231SHong Zhang for (j=0; j<anz; j++) { 1134a9d06231SHong Zhang row = adj[j]; 1135a9d06231SHong Zhang pnz = pi_loc[row+1] - pi_loc[row]; 1136a9d06231SHong Zhang pj = pj_loc + pi_loc[row]; 1137a9d06231SHong Zhang pa = pa_loc + pi_loc[row]; 1138a9d06231SHong Zhang 1139a9d06231SHong Zhang /* perform dense axpy */ 1140e900a757SHong Zhang valtmp = ada[j]; 1141a9d06231SHong Zhang for (k=0; k<pnz; k++) { 1142e900a757SHong Zhang apa[pj[k]] += valtmp*pa[k]; 1143a9d06231SHong Zhang } 1144a9d06231SHong Zhang ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr); 1145a9d06231SHong Zhang } 1146a9d06231SHong Zhang 1147a9d06231SHong Zhang /* off-diagonal portion of A */ 1148a9d06231SHong Zhang anz = aoi[i+1] - aoi[i]; 1149a9d06231SHong Zhang aoj = ao->j + aoi[i]; 1150a9d06231SHong Zhang aoa = ao->a + aoi[i]; 1151a9d06231SHong Zhang for (j=0; j<anz; j++) { 1152a9d06231SHong Zhang row = aoj[j]; 1153a9d06231SHong Zhang pnz = pi_oth[row+1] - pi_oth[row]; 1154a9d06231SHong Zhang pj = pj_oth + pi_oth[row]; 1155a9d06231SHong Zhang pa = pa_oth + pi_oth[row]; 1156a9d06231SHong Zhang 1157a9d06231SHong Zhang /* perform dense axpy */ 1158e900a757SHong Zhang valtmp = aoa[j]; 1159a9d06231SHong Zhang for (k=0; k<pnz; k++) { 1160e900a757SHong Zhang apa[pj[k]] += valtmp*pa[k]; 1161a9d06231SHong Zhang } 1162a9d06231SHong Zhang ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr); 1163a9d06231SHong Zhang } 1164b091e509SHong Zhang #if defined(PTAP_PROFILE) 11658563dfccSBarry Smith ierr = PetscTime(&t2_1);CHKERRQ(ierr); 1166b091e509SHong Zhang et2_AP += t2_1 - t2_0; 1167b091e509SHong Zhang #endif 1168a9d06231SHong Zhang 1169a9d06231SHong Zhang /* 2-b) Compute Cseq = P_loc[i,:]^T*AP[i,:] using outer product */ 1170a9d06231SHong Zhang /*--------------------------------------------------------------*/ 1171a9d06231SHong Zhang apnz = api[i+1] - api[i]; 1172a9d06231SHong Zhang /* put the value into Co=(p->B)^T*AP (off-diagonal part, send to others) */ 1173a9d06231SHong Zhang pnz = po->i[i+1] - po->i[i]; 1174e900a757SHong Zhang poJ = po->j + po->i[i]; 1175a9d06231SHong Zhang pA = po->a + po->i[i]; 1176a9d06231SHong Zhang for (j=0; j<pnz; j++) { 1177e900a757SHong Zhang row = poJ[j]; 1178e900a757SHong Zhang cnz = coi[row+1] - coi[row]; 1179e900a757SHong Zhang cj = coj + coi[row]; 1180e900a757SHong Zhang ca = coa + coi[row]; 1181a9d06231SHong Zhang /* perform dense axpy */ 1182e900a757SHong Zhang valtmp = pA[j]; 1183a9d06231SHong Zhang for (k=0; k<cnz; k++) { 1184e900a757SHong Zhang ca[k] += valtmp*apa[cj[k]]; 1185a9d06231SHong Zhang } 1186a9d06231SHong Zhang ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr); 1187a9d06231SHong Zhang } 118805d62848SHong Zhang #if 1 1189a9d06231SHong Zhang /* put the value into Cd (diagonal part) */ 1190a9d06231SHong Zhang pnz = pd->i[i+1] - pd->i[i]; 1191e900a757SHong Zhang pdJ = pd->j + pd->i[i]; 1192a9d06231SHong Zhang pA = pd->a + pd->i[i]; 1193a9d06231SHong Zhang for (j=0; j<pnz; j++) { 1194e900a757SHong Zhang row = pdJ[j]; 1195e900a757SHong Zhang cnz = bi[row+1] - bi[row]; 1196e900a757SHong Zhang cj = bj + bi[row]; 1197e900a757SHong Zhang ca = ba + bi[row]; 1198a9d06231SHong Zhang /* perform dense axpy */ 1199e900a757SHong Zhang valtmp = pA[j]; 1200a9d06231SHong Zhang for (k=0; k<cnz; k++) { 1201e900a757SHong Zhang ca[k] += valtmp*apa[cj[k]]; 1202a9d06231SHong Zhang } 1203a9d06231SHong Zhang ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr); 1204a9d06231SHong Zhang } 12059516a85cSHong Zhang #endif 1206a9d06231SHong Zhang /* zero the current row of A*P */ 1207a9d06231SHong Zhang for (k=0; k<apnz; k++) apa[apJ[k]] = 0.0; 1208b091e509SHong Zhang #if defined(PTAP_PROFILE) 12098563dfccSBarry Smith ierr = PetscTime(&t2_2);CHKERRQ(ierr); 121005d62848SHong Zhang ePtAP += t2_2 - t2_1; 1211b091e509SHong Zhang #endif 1212a9d06231SHong Zhang } 12139516a85cSHong Zhang 12149516a85cSHong Zhang if (rank == 100) { 12159516a85cSHong Zhang for (row=0; row<cm; row++) { 12169516a85cSHong Zhang printf("[%d] row %d: ",rank,row); 12179516a85cSHong Zhang cnz = bi[row+1] - bi[row]; 12189516a85cSHong Zhang for (j=0; j<cnz; j++) printf(" %g,",ba[bi[row]+j]); 12199516a85cSHong Zhang printf("\n"); 12209516a85cSHong Zhang } 12219516a85cSHong Zhang } 12229516a85cSHong Zhang 122338571addSHong Zhang } else { /* Do sparse axpy on apa (length of ap_rmax, stores A[i,:]*P) - scalable, but slower */ 1224b38d1a2dSJed Brown ierr = PetscInfo(C,"Using scalable sparse axpy\n");CHKERRQ(ierr); 122538571addSHong Zhang /*-----------------------------------------------------------------------------------------*/ 1226a9d06231SHong Zhang pA=pa_loc; 122742c57489SHong Zhang for (i=0; i<am; i++) { 1228b091e509SHong Zhang #if defined(PTAP_PROFILE) 12298563dfccSBarry Smith ierr = PetscTime(&t2_0);CHKERRQ(ierr); 1230b091e509SHong Zhang #endif 1231f4a743e1SHong Zhang /* form i-th sparse row of A*P */ 123282412ba7SHong Zhang apnz = api[i+1] - api[i]; 123382412ba7SHong Zhang apJ = apj + api[i]; 123442c57489SHong Zhang /* diagonal portion of A */ 123582412ba7SHong Zhang anz = adi[i+1] - adi[i]; 12361d633602SHong Zhang adj = ad->j + adi[i]; 12371d633602SHong Zhang ada = ad->a + adi[i]; 123882412ba7SHong Zhang for (j=0; j<anz; j++) { 12391d633602SHong Zhang row = adj[j]; 124082412ba7SHong Zhang pnz = pi_loc[row+1] - pi_loc[row]; 124182412ba7SHong Zhang pj = pj_loc + pi_loc[row]; 124282412ba7SHong Zhang pa = pa_loc + pi_loc[row]; 1243e900a757SHong Zhang valtmp = ada[j]; 1244f4a743e1SHong Zhang nextp = 0; 124582412ba7SHong Zhang for (k=0; nextp<pnz; k++) { 124682412ba7SHong Zhang if (apJ[k] == pj[nextp]) { /* col of AP == col of P */ 1247e900a757SHong Zhang apa[k] += valtmp*pa[nextp++]; 124842c57489SHong Zhang } 124942c57489SHong Zhang } 1250dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr); 125142c57489SHong Zhang } 125242c57489SHong Zhang /* off-diagonal portion of A */ 125382412ba7SHong Zhang anz = aoi[i+1] - aoi[i]; 12541d633602SHong Zhang aoj = ao->j + aoi[i]; 12551d633602SHong Zhang aoa = ao->a + aoi[i]; 125682412ba7SHong Zhang for (j=0; j<anz; j++) { 12571d633602SHong Zhang row = aoj[j]; 125882412ba7SHong Zhang pnz = pi_oth[row+1] - pi_oth[row]; 125982412ba7SHong Zhang pj = pj_oth + pi_oth[row]; 126082412ba7SHong Zhang pa = pa_oth + pi_oth[row]; 1261e900a757SHong Zhang valtmp = aoa[j]; 1262f4a743e1SHong Zhang nextp = 0; 126382412ba7SHong Zhang for (k=0; nextp<pnz; k++) { 126482412ba7SHong Zhang if (apJ[k] == pj[nextp]) { /* col of AP == col of P */ 1265e900a757SHong Zhang apa[k] += valtmp*pa[nextp++]; 126642c57489SHong Zhang } 126742c57489SHong Zhang } 1268dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr); 126942c57489SHong Zhang } 1270b091e509SHong Zhang #if defined(PTAP_PROFILE) 12718563dfccSBarry Smith ierr = PetscTime(&t2_1);CHKERRQ(ierr); 1272b091e509SHong Zhang et2_AP += t2_1 - t2_0; 1273b091e509SHong Zhang #endif 127442c57489SHong Zhang 1275a9d06231SHong Zhang /* 2-b) Compute Cseq = P_loc[i,:]^T*AP[i,:] using outer product */ 1276a9d06231SHong Zhang /*--------------------------------------------------------------*/ 127782412ba7SHong Zhang pnz = pi_loc[i+1] - pi_loc[i]; 1278f9473cf7SHong Zhang pJ = pj_loc + pi_loc[i]; 127982412ba7SHong Zhang for (j=0; j<pnz; j++) { 128042c57489SHong Zhang nextap = 0; 1281f9473cf7SHong Zhang row = pJ[j]; /* global index */ 128282412ba7SHong Zhang if (row < pcstart || row >=pcend) { /* put the value into Co */ 1283e900a757SHong Zhang row = *poJ; 1284e900a757SHong Zhang cj = coj + coi[row]; 1285e900a757SHong Zhang ca = coa + coi[row]; poJ++; 128682412ba7SHong Zhang } else { /* put the value into Cd */ 1287e900a757SHong Zhang row = *pdJ; 1288e900a757SHong Zhang cj = bj + bi[row]; 1289e900a757SHong Zhang ca = ba + bi[row]; pdJ++; 129042c57489SHong Zhang } 1291e900a757SHong Zhang valtmp = pA[j]; 129282412ba7SHong Zhang for (k=0; nextap<apnz; k++) { 1293e900a757SHong Zhang if (cj[k]==apJ[nextap]) ca[k] += valtmp*apa[nextap++]; 129442c57489SHong Zhang } 1295dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr); 1296de0260b3SHong Zhang } 12970496f32cSHong Zhang pA += pnz; 129842c57489SHong Zhang /* zero the current row info for A*P */ 129982412ba7SHong Zhang ierr = PetscMemzero(apa,apnz*sizeof(MatScalar));CHKERRQ(ierr); 1300b091e509SHong Zhang #if defined(PTAP_PROFILE) 13018563dfccSBarry Smith ierr = PetscTime(&t2_2);CHKERRQ(ierr); 130205d62848SHong Zhang ePtAP += t2_2 - t2_1; 1303b091e509SHong Zhang #endif 130442c57489SHong Zhang } 1305d9824c15SHong Zhang } 130638571addSHong Zhang #if defined(PTAP_PROFILE) 13078563dfccSBarry Smith ierr = PetscTime(&t2);CHKERRQ(ierr); 130838571addSHong Zhang #endif 130942c57489SHong Zhang 1310a9d06231SHong Zhang /* 3) send and recv matrix values coa */ 1311a9d06231SHong Zhang /*------------------------------------*/ 131242c57489SHong Zhang buf_ri = merge->buf_ri; 131342c57489SHong Zhang buf_rj = merge->buf_rj; 131442c57489SHong Zhang len_s = merge->len_s; 1315fc42d0c8SSatish Balay ierr = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr); 131642c57489SHong Zhang ierr = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr); 131742c57489SHong Zhang 1318dcca6d9dSJed Brown ierr = PetscMalloc2(merge->nsend+1,&s_waits,size,&status);CHKERRQ(ierr); 131942c57489SHong Zhang for (proc=0,k=0; proc<size; proc++) { 132042c57489SHong Zhang if (!len_s[proc]) continue; 1321de0260b3SHong Zhang i = merge->owners_co[proc]; 1322de0260b3SHong Zhang ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr); 132342c57489SHong Zhang k++; 132442c57489SHong Zhang } 13250c468ba9SBarry Smith if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);} 13260c468ba9SBarry Smith if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);} 132742c57489SHong Zhang 1328a9d06231SHong Zhang ierr = PetscFree2(s_waits,status);CHKERRQ(ierr); 132942c57489SHong Zhang ierr = PetscFree(r_waits);CHKERRQ(ierr); 133082412ba7SHong Zhang ierr = PetscFree(coa);CHKERRQ(ierr); 133138571addSHong Zhang #if defined(PTAP_PROFILE) 13328563dfccSBarry Smith ierr = PetscTime(&t3);CHKERRQ(ierr); 133338571addSHong Zhang #endif 133442c57489SHong Zhang 1335a9d06231SHong Zhang /* 4) insert local Cseq and received values into Cmpi */ 1336a9d06231SHong Zhang /*------------------------------------------------------*/ 1337dcca6d9dSJed Brown ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr); 133842c57489SHong Zhang for (k=0; k<merge->nrecv; k++) { 133942c57489SHong Zhang buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 134042c57489SHong Zhang nrows = *(buf_ri_k[k]); 134142c57489SHong Zhang nextrow[k] = buf_ri_k[k]+1; /* next row number of k-th recved i-structure */ 134282412ba7SHong Zhang nextci[k] = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure */ 134342c57489SHong Zhang } 134442c57489SHong Zhang 1345de0260b3SHong Zhang for (i=0; i<cm; i++) { 134682412ba7SHong Zhang row = owners[rank] + i; /* global row index of C_seq */ 134742c57489SHong Zhang bj_i = bj + bi[i]; /* col indices of the i-th row of C */ 1348de0260b3SHong Zhang ba_i = ba + bi[i]; 134982412ba7SHong Zhang bnz = bi[i+1] - bi[i]; 135042c57489SHong Zhang /* add received vals into ba */ 135142c57489SHong Zhang for (k=0; k<merge->nrecv; k++) { /* k-th received message */ 135242c57489SHong Zhang /* i-th row */ 135342c57489SHong Zhang if (i == *nextrow[k]) { 135482412ba7SHong Zhang cnz = *(nextci[k]+1) - *nextci[k]; 135582412ba7SHong Zhang cj = buf_rj[k] + *(nextci[k]); 135682412ba7SHong Zhang ca = abuf_r[k] + *(nextci[k]); 135782412ba7SHong Zhang nextcj = 0; 135882412ba7SHong Zhang for (j=0; nextcj<cnz; j++) { 135982412ba7SHong Zhang if (bj_i[j] == cj[nextcj]) { /* bcol == ccol */ 136082412ba7SHong Zhang ba_i[j] += ca[nextcj++]; 136142c57489SHong Zhang } 136242c57489SHong Zhang } 136382412ba7SHong Zhang nextrow[k]++; nextci[k]++; 1364c5af039cSHong Zhang ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr); 136542c57489SHong Zhang } 136642c57489SHong Zhang } 136782412ba7SHong Zhang ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr); 136842c57489SHong Zhang } 136942c57489SHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 137042c57489SHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 137142c57489SHong Zhang 1372de0260b3SHong Zhang ierr = PetscFree(ba);CHKERRQ(ierr); 1373c05d87d6SBarry Smith ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr); 137442c57489SHong Zhang ierr = PetscFree(abuf_r);CHKERRQ(ierr); 13750572522cSBarry Smith ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr); 137638571addSHong Zhang #if defined(PTAP_PROFILE) 13778563dfccSBarry Smith ierr = PetscTime(&t4);CHKERRQ(ierr); 13789516a85cSHong Zhang if (rank==1) { 137905d62848SHong Zhang ierr = PetscPrintf(MPI_COMM_SELF," [%d] PtAPNum %g/P + %g/PtAP( %g/A*P + %g/Pt*AP ) + %g/comm + %g/Cloc = %g\n\n",rank,eP,t2-t1,et2_AP,ePtAP,t3-t2,t4-t3,t4-t0);CHKERRQ(ierr); 13809516a85cSHong Zhang } 138138571addSHong Zhang #endif 138242c57489SHong Zhang PetscFunctionReturn(0); 138342c57489SHong Zhang } 1384