xref: /petsc/src/mat/impls/aij/mpi/mpiptap.c (revision 80bb463919f608f8e15c8f7bc4ef4eb3143dad43)
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);
350d3441aeSHong Zhang     ierr = MatDestroy(&ptap->AP_loc);CHKERRQ(ierr);
362259aa2eSHong Zhang     ierr = MatDestroy(&ptap->C_loc);CHKERRQ(ierr);
372259aa2eSHong Zhang     ierr = MatDestroy(&ptap->C_oth);CHKERRQ(ierr);
380d3441aeSHong Zhang 
39c5af039cSHong Zhang     if (ptap->api) {ierr = PetscFree(ptap->api);CHKERRQ(ierr);}
40c5af039cSHong Zhang     if (ptap->apj) {ierr = PetscFree(ptap->apj);CHKERRQ(ierr);}
41414894bdSHong Zhang     if (ptap->apa) {ierr = PetscFree(ptap->apa);CHKERRQ(ierr);}
42c8b0795cSMark F. Adams     if (merge) {
43cc6cb787SHong Zhang       ierr = PetscFree(merge->id_r);CHKERRQ(ierr);
44cc6cb787SHong Zhang       ierr = PetscFree(merge->len_s);CHKERRQ(ierr);
45cc6cb787SHong Zhang       ierr = PetscFree(merge->len_r);CHKERRQ(ierr);
46cc6cb787SHong Zhang       ierr = PetscFree(merge->bi);CHKERRQ(ierr);
47cc6cb787SHong Zhang       ierr = PetscFree(merge->bj);CHKERRQ(ierr);
48c05d87d6SBarry Smith       ierr = PetscFree(merge->buf_ri[0]);CHKERRQ(ierr);
49cc6cb787SHong Zhang       ierr = PetscFree(merge->buf_ri);CHKERRQ(ierr);
50c05d87d6SBarry Smith       ierr = PetscFree(merge->buf_rj[0]);CHKERRQ(ierr);
51cc6cb787SHong Zhang       ierr = PetscFree(merge->buf_rj);CHKERRQ(ierr);
5205b42c5fSBarry Smith       ierr = PetscFree(merge->coi);CHKERRQ(ierr);
5305b42c5fSBarry Smith       ierr = PetscFree(merge->coj);CHKERRQ(ierr);
5405b42c5fSBarry Smith       ierr = PetscFree(merge->owners_co);CHKERRQ(ierr);
556bf464f9SBarry Smith       ierr = PetscLayoutDestroy(&merge->rowmap);CHKERRQ(ierr);
56dce485f0SBarry Smith       ierr = merge->destroy(A);CHKERRQ(ierr);
57f8487c73SHong Zhang       ierr = PetscFree(ptap->merge);CHKERRQ(ierr);
58bf0cc555SLisandro Dalcin     }
59f8487c73SHong Zhang     ierr = PetscFree(ptap);CHKERRQ(ierr);
60c8b0795cSMark F. Adams   }
61cc6cb787SHong Zhang   PetscFunctionReturn(0);
62cc6cb787SHong Zhang }
63cc6cb787SHong Zhang 
64cc6cb787SHong Zhang #undef __FUNCT__
65cc6cb787SHong Zhang #define __FUNCT__ "MatDuplicate_MPIAIJ_MatPtAP"
66f0c0a3a6SBarry Smith PetscErrorCode MatDuplicate_MPIAIJ_MatPtAP(Mat A, MatDuplicateOption op, Mat *M)
67f0c0a3a6SBarry Smith {
68cc6cb787SHong Zhang   PetscErrorCode      ierr;
6911a6856fSHong Zhang   Mat_MPIAIJ          *a     = (Mat_MPIAIJ*)A->data;
7011a6856fSHong Zhang   Mat_PtAPMPI         *ptap  = a->ptap;
7111a6856fSHong Zhang   Mat_Merge_SeqsToMPI *merge = ptap->merge;
72cc6cb787SHong Zhang 
73cc6cb787SHong Zhang   PetscFunctionBegin;
74dce485f0SBarry Smith   ierr = (*merge->duplicate)(A,op,M);CHKERRQ(ierr);
752205254eSKarl Rupp 
7611a6856fSHong Zhang   (*M)->ops->destroy   = merge->destroy;
7711a6856fSHong Zhang   (*M)->ops->duplicate = merge->duplicate;
78cc6cb787SHong Zhang   PetscFunctionReturn(0);
79cc6cb787SHong Zhang }
80cc6cb787SHong Zhang 
8142c57489SHong Zhang #undef __FUNCT__
82cf3ca8ceSHong Zhang #define __FUNCT__ "MatPtAP_MPIAIJ_MPIAIJ"
83cf3ca8ceSHong Zhang PetscErrorCode MatPtAP_MPIAIJ_MPIAIJ(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C)
8442c57489SHong Zhang {
8542c57489SHong Zhang   PetscErrorCode ierr;
860d3441aeSHong Zhang   PetscBool      newalg=PETSC_FALSE;
8742c57489SHong Zhang 
8842c57489SHong Zhang   PetscFunctionBegin;
890d3441aeSHong Zhang   ierr = PetscOptionsGetBool(NULL,"-matptap_new",&newalg,NULL);CHKERRQ(ierr);
90cf3ca8ceSHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
913ff4c91cSHong Zhang     ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
920d3441aeSHong Zhang     if (newalg) {
930d3441aeSHong Zhang       ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ_new(A,P,fill,C);CHKERRQ(ierr);
940d3441aeSHong Zhang     } else {
95cf3ca8ceSHong Zhang       ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ(A,P,fill,C);CHKERRQ(ierr);
960d3441aeSHong Zhang     }
973ff4c91cSHong Zhang     ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
987a7894deSKris Buschelman   }
993ff4c91cSHong Zhang   ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
1000d3441aeSHong Zhang   if (newalg) {
1010d3441aeSHong Zhang     ierr = MatPtAPNumeric_MPIAIJ_MPIAIJ_new(A,P,*C);CHKERRQ(ierr);
1020d3441aeSHong Zhang   } else {
103cf3ca8ceSHong Zhang     ierr = MatPtAPNumeric_MPIAIJ_MPIAIJ(A,P,*C);CHKERRQ(ierr);
1040d3441aeSHong Zhang   }
1053ff4c91cSHong Zhang   ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
10642c57489SHong Zhang   PetscFunctionReturn(0);
10742c57489SHong Zhang }
10842c57489SHong Zhang 
10942c57489SHong Zhang #undef __FUNCT__
1100d3441aeSHong Zhang #define __FUNCT__ "MatPtAPSymbolic_MPIAIJ_MPIAIJ_new"
1110d3441aeSHong Zhang PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ_new(Mat A,Mat P,PetscReal fill,Mat *C)
1120d3441aeSHong Zhang {
1130d3441aeSHong Zhang   PetscErrorCode ierr;
1140d3441aeSHong Zhang   Mat_PtAPMPI    *ptap;
1152259aa2eSHong Zhang   Mat_MPIAIJ     *p=(Mat_MPIAIJ*)P->data;
1169ce11a7cSHong Zhang   Mat            AP;
1172259aa2eSHong Zhang   Mat_MPIAIJ     *c;
1182259aa2eSHong Zhang   MPI_Comm       comm;
1192259aa2eSHong Zhang   PetscMPIInt    size,rank;
120*80bb4639SHong Zhang   PetscLogDouble      t0,t1,t2,t3,t4;
121*80bb4639SHong Zhang 
122*80bb4639SHong Zhang   PetscFunctionBegin;
1232259aa2eSHong Zhang   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
1242259aa2eSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1252259aa2eSHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1262259aa2eSHong Zhang 
1272259aa2eSHong Zhang   /* check if matrix local sizes are compatible -- MV! */
1282259aa2eSHong Zhang   if (A->rmap->rstart != P->rmap->rstart || A->rmap->rend != P->rmap->rend) {
1292259aa2eSHong 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);
1302259aa2eSHong Zhang   }
1312259aa2eSHong Zhang   if (A->cmap->rstart != P->rmap->rstart || A->cmap->rend != P->rmap->rend) {
1322259aa2eSHong 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);
1332259aa2eSHong Zhang   }
1342259aa2eSHong Zhang 
13515a3b8e2SHong Zhang   //ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ(A,P,fill,C);CHKERRQ(ierr); //rm later !!!
13615a3b8e2SHong Zhang   //=============================================================================
13715a3b8e2SHong Zhang   Mat                 Cmpi;
13815a3b8e2SHong Zhang   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
139de817e96SHong Zhang   Mat_SeqAIJ          *p_loc;
140de817e96SHong Zhang   PetscInt            *pi_loc;
141de817e96SHong Zhang   PetscInt            *pdti,*pdtj,*poti,*potj,*ptJ,nnz;
14215a3b8e2SHong Zhang   PetscInt            *lnk,*owners_co,*coi,*coj,i,k,pnz,row;
14315a3b8e2SHong Zhang   PetscInt            am=A->rmap->n,pN=P->cmap->N,pm=P->rmap->n,pn=P->cmap->n;
14415a3b8e2SHong Zhang   PetscBT             lnkbt;
14515a3b8e2SHong Zhang   PetscMPIInt         tagi,tagj,*len_si,*len_s,*len_ri,icompleted=0;
14615a3b8e2SHong Zhang   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
14715a3b8e2SHong Zhang   PetscInt            len,proc,*dnz,*onz,*owners;
14815a3b8e2SHong Zhang   PetscInt            nzi,*pti,*ptj;
14915a3b8e2SHong Zhang   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
15015a3b8e2SHong Zhang   MPI_Request         *swaits,*rwaits;
15115a3b8e2SHong Zhang   MPI_Status          *sstatus,rstatus;
15215a3b8e2SHong Zhang   Mat_Merge_SeqsToMPI *merge;
15315a3b8e2SHong Zhang   PetscInt            *api,*apj,*Jptr,apnz,*prmap=p->garray,pon,nspacedouble=0,j,ap_rmax=0;
15415a3b8e2SHong Zhang   PetscReal           afill=1.0,afill_tmp;
15515a3b8e2SHong Zhang   PetscInt            rmax;
15615a3b8e2SHong Zhang 
157*80bb4639SHong Zhang   ierr = PetscTime(&t0);CHKERRQ(ierr);
15815a3b8e2SHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
15915a3b8e2SHong Zhang   ierr        = PetscNew(&ptap);CHKERRQ(ierr);
16015a3b8e2SHong Zhang   ierr        = PetscNew(&merge);CHKERRQ(ierr);
16115a3b8e2SHong Zhang   ptap->merge = merge;
16215a3b8e2SHong Zhang   ptap->reuse = MAT_INITIAL_MATRIX;
16315a3b8e2SHong Zhang 
16415a3b8e2SHong Zhang   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
16515a3b8e2SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
16615a3b8e2SHong Zhang 
16715a3b8e2SHong Zhang   /* get P_loc by taking all local rows of P */
16815a3b8e2SHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
16915a3b8e2SHong Zhang 
170de817e96SHong Zhang   ptap->reuse = MAT_INITIAL_MATRIX;
17115a3b8e2SHong Zhang 
172de817e96SHong Zhang   /* (1) compute symbolic AP = A*P, then get AP_loc */
17315a3b8e2SHong Zhang   /*--------------------------------------------------------------------------*/
174de817e96SHong Zhang   ierr = MatTranspose_SeqAIJ(p->A,MAT_INITIAL_MATRIX,&ptap->Rd);CHKERRQ(ierr);
175de817e96SHong Zhang   ierr = MatTranspose_SeqAIJ(p->B,MAT_INITIAL_MATRIX,&ptap->Ro);CHKERRQ(ierr);
17615a3b8e2SHong Zhang 
177de817e96SHong Zhang   ierr = MatMatMult(A,P,MAT_INITIAL_MATRIX,2.0,&AP);CHKERRQ(ierr);
178de817e96SHong Zhang   ierr = MatMPIAIJGetLocalMat(AP,MAT_INITIAL_MATRIX,&ptap->AP_loc);CHKERRQ(ierr);
179de817e96SHong Zhang   ierr = MatDestroy(&AP);CHKERRQ(ierr);
18015a3b8e2SHong Zhang 
181de817e96SHong Zhang   /* (2) compute C_loc=Rd*AP_loc, Co=Ro*AP_loc */
182de817e96SHong Zhang   ierr = MatMatMult_SeqAIJ_SeqAIJ(ptap->Rd,ptap->AP_loc,MAT_INITIAL_MATRIX,2.0,&ptap->C_loc);CHKERRQ(ierr);
183de817e96SHong Zhang   ierr = MatMatMult_SeqAIJ_SeqAIJ(ptap->Ro,ptap->AP_loc,MAT_INITIAL_MATRIX,2.0,&ptap->C_oth);CHKERRQ(ierr);
184*80bb4639SHong Zhang   ierr = PetscTime(&t1);CHKERRQ(ierr);
18515a3b8e2SHong Zhang 
186de817e96SHong Zhang   p_loc  = (Mat_SeqAIJ*)(ptap->P_loc)->data;
187de817e96SHong Zhang   pi_loc = p_loc->i;
18815a3b8e2SHong Zhang 
189de817e96SHong Zhang   Mat_SeqAIJ *ap = (Mat_SeqAIJ*)ptap->AP_loc->data;
190de817e96SHong Zhang   api = ap->i;
191de817e96SHong Zhang   apj = ap->j;
19215a3b8e2SHong Zhang 
19315a3b8e2SHong Zhang   /* (2) determine symbolic Co=(p->B)^T*AP - send to others (coi,coj)*/
19415a3b8e2SHong Zhang   /*-----------------------------------------------------------------*/
195*80bb4639SHong Zhang   Mat_SeqAIJ *po = (Mat_SeqAIJ *)ptap->Ro->data;
196*80bb4639SHong Zhang   poti = po->i; potj = po->j;
19715a3b8e2SHong Zhang 
19815a3b8e2SHong Zhang   /* then, compute symbolic Co = (p->B)^T*AP */
19915a3b8e2SHong Zhang   pon    = (p->B)->cmap->n; /* total num of rows to be sent to other processors
20015a3b8e2SHong Zhang                          >= (num of nonzero rows of C_seq) - pn */
20115a3b8e2SHong Zhang   ierr   = PetscMalloc1(pon+1,&coi);CHKERRQ(ierr);
20215a3b8e2SHong Zhang   coi[0] = 0;
20315a3b8e2SHong Zhang 
20415a3b8e2SHong Zhang   /* set initial free space to be fill*(nnz(p->B) + nnz(AP)) */
20515a3b8e2SHong Zhang   nnz           = fill*(poti[pon] + api[am]);
20615a3b8e2SHong Zhang   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
20715a3b8e2SHong Zhang   current_space = free_space;
20815a3b8e2SHong Zhang 
209de817e96SHong Zhang   /* create and initialize a linked list */
210de817e96SHong Zhang   ierr = PetscLLCondensedCreate(pN,pN,&lnk,&lnkbt);CHKERRQ(ierr);
211de817e96SHong Zhang 
21215a3b8e2SHong Zhang   for (i=0; i<pon; i++) {
21315a3b8e2SHong Zhang     pnz = poti[i+1] - poti[i];
21415a3b8e2SHong Zhang     ptJ = potj + poti[i];
21515a3b8e2SHong Zhang     for (j=0; j<pnz; j++) {
21615a3b8e2SHong Zhang       row  = ptJ[j]; /* row of AP == col of Pot */
21715a3b8e2SHong Zhang       apnz = api[row+1] - api[row];
21815a3b8e2SHong Zhang       Jptr = apj + api[row];
21915a3b8e2SHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
22015a3b8e2SHong Zhang       ierr = PetscLLCondensedAddSorted(apnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
22115a3b8e2SHong Zhang     }
22215a3b8e2SHong Zhang     nnz = lnk[0];
22315a3b8e2SHong Zhang 
22415a3b8e2SHong Zhang     /* If free space is not available, double the total space in the list */
22515a3b8e2SHong Zhang     if (current_space->local_remaining<nnz) {
22615a3b8e2SHong Zhang       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
22715a3b8e2SHong Zhang       nspacedouble++;
22815a3b8e2SHong Zhang     }
22915a3b8e2SHong Zhang 
23015a3b8e2SHong Zhang     /* Copy data into free space, and zero out denserows */
23115a3b8e2SHong Zhang     ierr = PetscLLCondensedClean(pN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
23215a3b8e2SHong Zhang 
23315a3b8e2SHong Zhang     current_space->array           += nnz;
23415a3b8e2SHong Zhang     current_space->local_used      += nnz;
23515a3b8e2SHong Zhang     current_space->local_remaining -= nnz;
23615a3b8e2SHong Zhang 
23715a3b8e2SHong Zhang     coi[i+1] = coi[i] + nnz;
23815a3b8e2SHong Zhang   }
23915a3b8e2SHong Zhang 
24015a3b8e2SHong Zhang   ierr      = PetscMalloc1(coi[pon],&coj);CHKERRQ(ierr);
24115a3b8e2SHong Zhang   ierr      = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr);
24215a3b8e2SHong Zhang   afill_tmp = (PetscReal)coi[pon]/(poti[pon] + api[am]+1);
24315a3b8e2SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
24415a3b8e2SHong Zhang 
24515a3b8e2SHong Zhang   /* (3) send j-array (coj) of Co to other processors */
24615a3b8e2SHong Zhang   /*--------------------------------------------------*/
24715a3b8e2SHong Zhang   ierr = PetscCalloc1(size,&merge->len_s);CHKERRQ(ierr);
24815a3b8e2SHong Zhang   len_s        = merge->len_s;
24915a3b8e2SHong Zhang   merge->nsend = 0;
25015a3b8e2SHong Zhang 
25115a3b8e2SHong Zhang   /* determine row ownership */
25215a3b8e2SHong Zhang   ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr);
25315a3b8e2SHong Zhang   merge->rowmap->n  = pn;
25415a3b8e2SHong Zhang   merge->rowmap->bs = 1;
25515a3b8e2SHong Zhang 
25615a3b8e2SHong Zhang   ierr   = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr);
25715a3b8e2SHong Zhang   owners = merge->rowmap->range;
25815a3b8e2SHong Zhang 
25915a3b8e2SHong Zhang   /* determine the number of messages to send, their lengths */
26015a3b8e2SHong Zhang   ierr = PetscMalloc2(size,&len_si,size,&sstatus);CHKERRQ(ierr);
26115a3b8e2SHong Zhang   ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
26215a3b8e2SHong Zhang   ierr = PetscMalloc1(size+2,&owners_co);CHKERRQ(ierr);
26315a3b8e2SHong Zhang 
26415a3b8e2SHong Zhang   proc = 0;
26515a3b8e2SHong Zhang   for (i=0; i<pon; i++) {
26615a3b8e2SHong Zhang     while (prmap[i] >= owners[proc+1]) proc++;
26715a3b8e2SHong Zhang     len_si[proc]++;               /* num of rows in Co(=Pt*AP) to be sent to [proc] */
26815a3b8e2SHong Zhang     len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */
26915a3b8e2SHong Zhang   }
27015a3b8e2SHong Zhang 
27115a3b8e2SHong Zhang   len          = 0; /* max length of buf_si[], see (4) */
27215a3b8e2SHong Zhang   owners_co[0] = 0;
27315a3b8e2SHong Zhang   for (proc=0; proc<size; proc++) {
27415a3b8e2SHong Zhang     owners_co[proc+1] = owners_co[proc] + len_si[proc];
27515a3b8e2SHong Zhang     if (len_s[proc]) {
27615a3b8e2SHong Zhang       merge->nsend++;
27715a3b8e2SHong Zhang       len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */
27815a3b8e2SHong Zhang       len         += len_si[proc];
27915a3b8e2SHong Zhang     }
28015a3b8e2SHong Zhang   }
28115a3b8e2SHong Zhang 
28215a3b8e2SHong Zhang   /* determine the number and length of messages to receive for coi and coj  */
28315a3b8e2SHong Zhang   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&merge->nrecv);CHKERRQ(ierr);
28415a3b8e2SHong Zhang   ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr);
28515a3b8e2SHong Zhang 
28615a3b8e2SHong Zhang   /* post the Irecv and Isend of coj */
28715a3b8e2SHong Zhang   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
28815a3b8e2SHong Zhang   ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
28915a3b8e2SHong Zhang   ierr = PetscMalloc1(merge->nsend+1,&swaits);CHKERRQ(ierr);
29015a3b8e2SHong Zhang   for (proc=0, k=0; proc<size; proc++) {
29115a3b8e2SHong Zhang     if (!len_s[proc]) continue;
29215a3b8e2SHong Zhang     i    = owners_co[proc];
29315a3b8e2SHong Zhang     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
29415a3b8e2SHong Zhang     k++;
29515a3b8e2SHong Zhang   }
29615a3b8e2SHong Zhang 
29715a3b8e2SHong Zhang   /* receives and sends of coj are complete */
29815a3b8e2SHong Zhang   for (i=0; i<merge->nrecv; i++) {
29915a3b8e2SHong Zhang     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
30015a3b8e2SHong Zhang   }
30115a3b8e2SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
30215a3b8e2SHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
30315a3b8e2SHong Zhang 
30415a3b8e2SHong Zhang   /* (4) send and recv coi */
30515a3b8e2SHong Zhang   /*-----------------------*/
30615a3b8e2SHong Zhang   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
30715a3b8e2SHong Zhang   ierr   = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
30815a3b8e2SHong Zhang   ierr   = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr);
30915a3b8e2SHong Zhang   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
31015a3b8e2SHong Zhang   for (proc=0,k=0; proc<size; proc++) {
31115a3b8e2SHong Zhang     if (!len_s[proc]) continue;
31215a3b8e2SHong Zhang     /* form outgoing message for i-structure:
31315a3b8e2SHong Zhang          buf_si[0]:                 nrows to be sent
31415a3b8e2SHong Zhang                [1:nrows]:           row index (global)
31515a3b8e2SHong Zhang                [nrows+1:2*nrows+1]: i-structure index
31615a3b8e2SHong Zhang     */
31715a3b8e2SHong Zhang     /*-------------------------------------------*/
31815a3b8e2SHong Zhang     nrows       = len_si[proc]/2 - 1; /* num of rows in Co to be sent to [proc] */
31915a3b8e2SHong Zhang     buf_si_i    = buf_si + nrows+1;
32015a3b8e2SHong Zhang     buf_si[0]   = nrows;
32115a3b8e2SHong Zhang     buf_si_i[0] = 0;
32215a3b8e2SHong Zhang     nrows       = 0;
32315a3b8e2SHong Zhang     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
32415a3b8e2SHong Zhang       nzi = coi[i+1] - coi[i];
32515a3b8e2SHong Zhang       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi;  /* i-structure */
32615a3b8e2SHong Zhang       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
32715a3b8e2SHong Zhang       nrows++;
32815a3b8e2SHong Zhang     }
32915a3b8e2SHong Zhang     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
33015a3b8e2SHong Zhang     k++;
33115a3b8e2SHong Zhang     buf_si += len_si[proc];
33215a3b8e2SHong Zhang   }
33315a3b8e2SHong Zhang   i = merge->nrecv;
33415a3b8e2SHong Zhang   while (i--) {
33515a3b8e2SHong Zhang     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
33615a3b8e2SHong Zhang   }
33715a3b8e2SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
33815a3b8e2SHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
33915a3b8e2SHong Zhang 
34015a3b8e2SHong Zhang   ierr = PetscFree2(len_si,sstatus);CHKERRQ(ierr);
34115a3b8e2SHong Zhang   ierr = PetscFree(len_ri);CHKERRQ(ierr);
34215a3b8e2SHong Zhang   ierr = PetscFree(swaits);CHKERRQ(ierr);
34315a3b8e2SHong Zhang   ierr = PetscFree(buf_s);CHKERRQ(ierr);
34415a3b8e2SHong Zhang 
345*80bb4639SHong Zhang   ierr = PetscTime(&t2);CHKERRQ(ierr);
346*80bb4639SHong Zhang 
34715a3b8e2SHong Zhang   /* (5) compute the local portion of C (mpi mat) */
34815a3b8e2SHong Zhang   /*----------------------------------------------*/
349*80bb4639SHong Zhang   Mat_SeqAIJ *pd = (Mat_SeqAIJ *)ptap->Rd->data;
350*80bb4639SHong Zhang   pdti = pd->i; pdtj = pd->j;
35115a3b8e2SHong Zhang 
35215a3b8e2SHong Zhang   /* allocate pti array and free space for accumulating nonzero column info */
35315a3b8e2SHong Zhang   ierr   = PetscMalloc1(pn+1,&pti);CHKERRQ(ierr);
35415a3b8e2SHong Zhang   pti[0] = 0;
35515a3b8e2SHong Zhang 
35615a3b8e2SHong Zhang   /* set initial free space to be fill*(nnz(P) + nnz(AP)) */
35715a3b8e2SHong Zhang   nnz           = fill*(pi_loc[pm] + api[am]);
35815a3b8e2SHong Zhang   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
35915a3b8e2SHong Zhang   current_space = free_space;
36015a3b8e2SHong Zhang 
36115a3b8e2SHong Zhang   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
36215a3b8e2SHong Zhang   for (k=0; k<merge->nrecv; k++) {
36315a3b8e2SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
36415a3b8e2SHong Zhang     nrows       = *buf_ri_k[k];
36515a3b8e2SHong Zhang     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
36615a3b8e2SHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
36715a3b8e2SHong Zhang   }
36815a3b8e2SHong Zhang   ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr);
36915a3b8e2SHong Zhang   rmax = 0;
37015a3b8e2SHong Zhang   for (i=0; i<pn; i++) {
37115a3b8e2SHong Zhang     /* add pdt[i,:]*AP into lnk */
37215a3b8e2SHong Zhang     pnz = pdti[i+1] - pdti[i];
37315a3b8e2SHong Zhang     ptJ = pdtj + pdti[i];
37415a3b8e2SHong Zhang     for (j=0; j<pnz; j++) {
37515a3b8e2SHong Zhang       row  = ptJ[j];  /* row of AP == col of Pt */
37615a3b8e2SHong Zhang       apnz = api[row+1] - api[row];
37715a3b8e2SHong Zhang       Jptr = apj + api[row];
37815a3b8e2SHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
37915a3b8e2SHong Zhang       ierr = PetscLLCondensedAddSorted(apnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
38015a3b8e2SHong Zhang     }
38115a3b8e2SHong Zhang 
38215a3b8e2SHong Zhang     /* add received col data into lnk */
38315a3b8e2SHong Zhang     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
38415a3b8e2SHong Zhang       if (i == *nextrow[k]) { /* i-th row */
38515a3b8e2SHong Zhang         nzi  = *(nextci[k]+1) - *nextci[k];
38615a3b8e2SHong Zhang         Jptr = buf_rj[k] + *nextci[k];
38715a3b8e2SHong Zhang         ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr);
38815a3b8e2SHong Zhang         nextrow[k]++; nextci[k]++;
38915a3b8e2SHong Zhang       }
39015a3b8e2SHong Zhang     }
39115a3b8e2SHong Zhang     nnz = lnk[0];
39215a3b8e2SHong Zhang 
39315a3b8e2SHong Zhang     /* if free space is not available, make more free space */
39415a3b8e2SHong Zhang     if (current_space->local_remaining<nnz) {
39515a3b8e2SHong Zhang       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
39615a3b8e2SHong Zhang       nspacedouble++;
39715a3b8e2SHong Zhang     }
39815a3b8e2SHong Zhang     /* copy data into free space, then initialize lnk */
39915a3b8e2SHong Zhang     ierr = PetscLLCondensedClean(pN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
40015a3b8e2SHong Zhang     ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr);
40115a3b8e2SHong Zhang 
40215a3b8e2SHong Zhang     current_space->array           += nnz;
40315a3b8e2SHong Zhang     current_space->local_used      += nnz;
40415a3b8e2SHong Zhang     current_space->local_remaining -= nnz;
40515a3b8e2SHong Zhang 
40615a3b8e2SHong Zhang     pti[i+1] = pti[i] + nnz;
40715a3b8e2SHong Zhang     if (nnz > rmax) rmax = nnz;
40815a3b8e2SHong Zhang   }
40915a3b8e2SHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
41015a3b8e2SHong Zhang 
41115a3b8e2SHong Zhang   ierr      = PetscMalloc1(pti[pn]+1,&ptj);CHKERRQ(ierr);
41215a3b8e2SHong Zhang   ierr      = PetscFreeSpaceContiguous(&free_space,ptj);CHKERRQ(ierr);
41315a3b8e2SHong Zhang   afill_tmp = (PetscReal)pti[pn]/(pi_loc[pm] + api[am]+1);
41415a3b8e2SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
41515a3b8e2SHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
41615a3b8e2SHong Zhang 
417*80bb4639SHong Zhang   ierr = PetscTime(&t3);CHKERRQ(ierr);
418*80bb4639SHong Zhang 
41915a3b8e2SHong Zhang   /* (6) create symbolic parallel matrix Cmpi */
42015a3b8e2SHong Zhang   /*------------------------------------------*/
42115a3b8e2SHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
42215a3b8e2SHong Zhang   ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
42315a3b8e2SHong Zhang   ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr);
42415a3b8e2SHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
42515a3b8e2SHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
42615a3b8e2SHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
42715a3b8e2SHong Zhang 
42815a3b8e2SHong Zhang   merge->bi        = pti;      /* Cseq->i */
42915a3b8e2SHong Zhang   merge->bj        = ptj;      /* Cseq->j */
43015a3b8e2SHong Zhang   merge->coi       = coi;      /* Co->i   */
43115a3b8e2SHong Zhang   merge->coj       = coj;      /* Co->j   */
43215a3b8e2SHong Zhang   merge->buf_ri    = buf_ri;
43315a3b8e2SHong Zhang   merge->buf_rj    = buf_rj;
43415a3b8e2SHong Zhang   merge->owners_co = owners_co;
43515a3b8e2SHong Zhang   merge->destroy   = Cmpi->ops->destroy;
43615a3b8e2SHong Zhang   merge->duplicate = Cmpi->ops->duplicate;
43715a3b8e2SHong Zhang 
43815a3b8e2SHong Zhang   /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */
43915a3b8e2SHong Zhang   Cmpi->assembled      = PETSC_FALSE;
44015a3b8e2SHong Zhang   Cmpi->ops->destroy   = MatDestroy_MPIAIJ_PtAP;
44115a3b8e2SHong Zhang   Cmpi->ops->duplicate = MatDuplicate_MPIAIJ_MatPtAP;
44215a3b8e2SHong Zhang 
44315a3b8e2SHong Zhang   /* attach the supporting struct to Cmpi for reuse */
44415a3b8e2SHong Zhang   c           = (Mat_MPIAIJ*)Cmpi->data;
44515a3b8e2SHong Zhang   c->ptap     = ptap;
446de817e96SHong Zhang   ptap->api   = NULL;
447de817e96SHong Zhang   ptap->apj   = NULL;
44815a3b8e2SHong Zhang   ptap->rmax  = ap_rmax;
44915a3b8e2SHong Zhang   *C          = Cmpi;
45015a3b8e2SHong Zhang 
4512259aa2eSHong Zhang   /* flag 'scalable' determines which implementations to be used:
4522259aa2eSHong Zhang        0: do dense axpy in MatPtAPNumeric() - fast, but requires storage of a nonscalable dense array apa;
4532259aa2eSHong Zhang        1: do sparse axpy in MatPtAPNumeric() - might slow, uses a sparse array apa */
4542259aa2eSHong Zhang   /* set default scalable */
4552259aa2eSHong Zhang   ptap->scalable = PETSC_FALSE; //PETSC_TRUE;
4562259aa2eSHong Zhang 
45715a3b8e2SHong Zhang   ierr = PetscOptionsGetBool(((PetscObject)(*C))->prefix,"-matptap_scalable",&ptap->scalable,NULL);CHKERRQ(ierr);
4582259aa2eSHong Zhang   if (!ptap->scalable) {  /* Do dense axpy */
4592259aa2eSHong Zhang     ierr = PetscCalloc1(P->cmap->N,&ptap->apa);CHKERRQ(ierr);
4602259aa2eSHong Zhang   } else {
4612259aa2eSHong Zhang     //ierr = PetscCalloc1(ap_rmax+1,&ptap->apa);CHKERRQ(ierr);
4622259aa2eSHong Zhang   }
4632259aa2eSHong Zhang 
464*80bb4639SHong Zhang   ierr = PetscTime(&t4);CHKERRQ(ierr);
465*80bb4639SHong Zhang   if (rank == 1) printf("PtAPSym: %g + %g + %g = %g\n",t1-t0,t2-t1,t3-t2,t4-t3);
466*80bb4639SHong Zhang 
4670d3441aeSHong Zhang   PetscFunctionReturn(0);
4680d3441aeSHong Zhang }
4690d3441aeSHong Zhang 
4700d3441aeSHong Zhang #undef __FUNCT__
4710d3441aeSHong Zhang #define __FUNCT__ "MatPtAPNumeric_MPIAIJ_MPIAIJ_new"
4720d3441aeSHong Zhang PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ_new(Mat A,Mat P,Mat C)
4730d3441aeSHong Zhang {
4740d3441aeSHong Zhang   PetscErrorCode    ierr;
4750d3441aeSHong Zhang   Mat_MPIAIJ        *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
4760d3441aeSHong Zhang   Mat_SeqAIJ        *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
4770d3441aeSHong Zhang   Mat_PtAPMPI       *ptap = c->ptap;
4789ce11a7cSHong Zhang   Mat               AP_loc,C_loc,C_oth;
4790d3441aeSHong Zhang   PetscInt          i,rstart,rend,cm,ncols,row;
4800d3441aeSHong Zhang   PetscMPIInt       rank;
4810d3441aeSHong Zhang   MPI_Comm          comm;
4820d3441aeSHong Zhang   const PetscInt    *cols;
4830d3441aeSHong Zhang   const PetscScalar *vals;
4840d3441aeSHong Zhang   PetscLogDouble    t0,t1,t2,t3,t4,eR,eAP,eCseq,eCmpi;
4850d3441aeSHong Zhang 
4860d3441aeSHong Zhang   PetscFunctionBegin;
4870d3441aeSHong Zhang   ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr);
4880d3441aeSHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
4890d3441aeSHong Zhang 
4900d3441aeSHong Zhang   ierr = MatZeroEntries(C);CHKERRQ(ierr);
4910d3441aeSHong Zhang 
492e497d3c8SHong Zhang   /* 1) get R = Pd^T,Ro = Po^T */
4930d3441aeSHong Zhang   ierr = PetscTime(&t0);CHKERRQ(ierr);
4940d3441aeSHong Zhang   ierr = MatTranspose_SeqAIJ(p->A,MAT_REUSE_MATRIX,&ptap->Rd);CHKERRQ(ierr);
4950d3441aeSHong Zhang   ierr = MatTranspose_SeqAIJ(p->B,MAT_REUSE_MATRIX,&ptap->Ro);CHKERRQ(ierr);
4960d3441aeSHong Zhang   ierr = PetscTime(&t1);CHKERRQ(ierr);
4970d3441aeSHong Zhang   eR = t1 - t0;
4980d3441aeSHong Zhang 
499e497d3c8SHong Zhang   /* 2) get AP_loc */
5000d3441aeSHong Zhang   AP_loc = ptap->AP_loc;
501e497d3c8SHong Zhang   Mat_SeqAIJ    *ap=(Mat_SeqAIJ*)AP_loc->data;
5020d3441aeSHong Zhang 
503e497d3c8SHong Zhang   /* 2-1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
5040d3441aeSHong Zhang   /*-----------------------------------------------------*/
505e497d3c8SHong Zhang   if (ptap->reuse == MAT_INITIAL_MATRIX) {
5060d3441aeSHong Zhang     /* P_oth and P_loc are obtained in MatPtASymbolic(), skip calling MatGetBrowsOfAoCols() and MatMPIAIJGetLocalMat() */
507e497d3c8SHong Zhang     ptap->reuse = MAT_REUSE_MATRIX;
508e497d3c8SHong Zhang   } else { /* update numerical values of P_oth and P_loc */
5090d3441aeSHong Zhang     ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
5100d3441aeSHong Zhang     ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
511e497d3c8SHong Zhang   }
5120d3441aeSHong Zhang 
513e497d3c8SHong Zhang   /* 2-2) compute numeric C_seq = P_loc^T*A_loc*P - dominating part */
5140d3441aeSHong Zhang   /*--------------------------------------------------------------*/
5150d3441aeSHong Zhang   /* get data from symbolic products */
5160d3441aeSHong Zhang   Mat_SeqAIJ  *p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
5170d3441aeSHong Zhang   Mat_SeqAIJ  *p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
518e497d3c8SHong Zhang   PetscInt    *api,*apj,am = A->rmap->n,j,col,apnz;
519e497d3c8SHong Zhang   PetscScalar *apa = ptap->apa;
5200d3441aeSHong Zhang 
5212259aa2eSHong Zhang   api = ap->i;
5222259aa2eSHong Zhang   apj = ap->j;
523e497d3c8SHong Zhang   for (i=0; i<am; i++) {
524e497d3c8SHong Zhang     /* AP[i,:] = A[i,:]*P = Ad*P_loc Ao*P_oth */
525e497d3c8SHong Zhang     AProw_nonscalable(i,ad,ao,p_loc,p_oth,apa);
526e497d3c8SHong Zhang     apnz = api[i+1] - api[i];
527e497d3c8SHong Zhang     for (j=0; j<apnz; j++) {
528e497d3c8SHong Zhang       col = apj[j+api[i]];
529e497d3c8SHong Zhang       ap->a[j+ap->i[i]] = apa[col];
530e497d3c8SHong Zhang       apa[col] = 0.0;
5310d3441aeSHong Zhang     }
5320d3441aeSHong Zhang   }
5330d3441aeSHong Zhang 
5340d3441aeSHong Zhang   ierr = PetscTime(&t2);CHKERRQ(ierr);
5350d3441aeSHong Zhang   eAP = t2 - t1;
5360d3441aeSHong Zhang 
537e497d3c8SHong Zhang   /* 3) C_loc = R*AP_loc, Co = Ro*AP_loc */
5389ce11a7cSHong Zhang   ierr = MatMatMult_SeqAIJ_SeqAIJ(ptap->Rd,AP_loc,MAT_REUSE_MATRIX,2.0,&ptap->C_loc);CHKERRQ(ierr);
5399ce11a7cSHong Zhang   ierr = MatMatMult_SeqAIJ_SeqAIJ(ptap->Ro,AP_loc,MAT_REUSE_MATRIX,2.0,&ptap->C_oth);CHKERRQ(ierr);
5409ce11a7cSHong Zhang   C_loc = ptap->C_loc;
5419ce11a7cSHong Zhang   C_oth = ptap->C_oth;
5420d3441aeSHong Zhang   //printf("[%d] Co %d, %d\n", rank,Co->rmap->N,Co->cmap->N);
5430d3441aeSHong Zhang   ierr = PetscTime(&t3);CHKERRQ(ierr);
5440d3441aeSHong Zhang   eCseq = t3 - t2;
5450d3441aeSHong Zhang 
5460d3441aeSHong Zhang   /* add C_loc and Co to to C */
5470d3441aeSHong Zhang   ierr = MatGetOwnershipRange(C,&rstart,&rend);CHKERRQ(ierr);
5480d3441aeSHong Zhang 
5490d3441aeSHong Zhang   /* C_loc -> C */
5500d3441aeSHong Zhang   cm = C_loc->rmap->N;
5519ce11a7cSHong Zhang   Mat_SeqAIJ *c_seq;
5529ce11a7cSHong Zhang   c_seq = (Mat_SeqAIJ*)C_loc->data;
5530d3441aeSHong Zhang   for (i=0; i<cm; i++) {
5549ce11a7cSHong Zhang     ncols = c_seq->i[i+1] - c_seq->i[i];
5550d3441aeSHong Zhang     row = rstart + i;
5569ce11a7cSHong Zhang     cols = c_seq->j + c_seq->i[i];
5579ce11a7cSHong Zhang     vals = c_seq->a + c_seq->i[i];
5580d3441aeSHong Zhang     ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr);
5590d3441aeSHong Zhang   }
5600d3441aeSHong Zhang 
5610d3441aeSHong Zhang   /* Co -> C, off-processor part */
5620d3441aeSHong Zhang   //printf("[%d] p->B %d, %d\n",rank,p->B->rmap->N,p->B->cmap->N);
5639ce11a7cSHong Zhang   cm = C_oth->rmap->N;
5649ce11a7cSHong Zhang   c_seq = (Mat_SeqAIJ*)C_oth->data;
5659ce11a7cSHong Zhang   for (i=0; i<cm; i++) {
5669ce11a7cSHong Zhang     ncols = c_seq->i[i+1] - c_seq->i[i];
5670d3441aeSHong Zhang     row = p->garray[i];
5689ce11a7cSHong Zhang     cols = c_seq->j + c_seq->i[i];
5699ce11a7cSHong Zhang     vals = c_seq->a + c_seq->i[i];
5700d3441aeSHong Zhang     //printf("[%d] row[%d] = %d\n",rank,i,row);
5710d3441aeSHong Zhang     ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr);
5720d3441aeSHong Zhang   }
5730d3441aeSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
5740d3441aeSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
5750d3441aeSHong Zhang   ierr = PetscTime(&t4);CHKERRQ(ierr);
5760d3441aeSHong Zhang   eCmpi = t4 - t3;
5770d3441aeSHong Zhang 
5780d3441aeSHong Zhang   if (rank==1) {
5790d3441aeSHong 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);
5800d3441aeSHong Zhang   }
5810d3441aeSHong Zhang   PetscFunctionReturn(0);
5820d3441aeSHong Zhang }
5830d3441aeSHong Zhang 
5840d3441aeSHong Zhang #undef __FUNCT__
58542c57489SHong Zhang #define __FUNCT__ "MatPtAPSymbolic_MPIAIJ_MPIAIJ"
58642c57489SHong Zhang PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
58742c57489SHong Zhang {
58842c57489SHong Zhang   PetscErrorCode      ierr;
58977584fe6SHong Zhang   Mat                 Cmpi;
590f8487c73SHong Zhang   Mat_PtAPMPI         *ptap;
5910298fd71SBarry Smith   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
592f8487c73SHong Zhang   Mat_MPIAIJ          *a        =(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c;
59342c57489SHong Zhang   Mat_SeqAIJ          *ad       =(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
59442c57489SHong Zhang   Mat_SeqAIJ          *p_loc,*p_oth;
595ded4f561SHong Zhang   PetscInt            *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pdti,*pdtj,*poti,*potj,*ptJ;
59677584fe6SHong Zhang   PetscInt            *adi=ad->i,*aj,*aoi=ao->i,nnz;
597a914927fSHong Zhang   PetscInt            *lnk,*owners_co,*coi,*coj,i,k,pnz,row;
598d16ca5f0SHong Zhang   PetscInt            am=A->rmap->n,pN=P->cmap->N,pm=P->rmap->n,pn=P->cmap->n;
59942c57489SHong Zhang   PetscBT             lnkbt;
600ce94432eSBarry Smith   MPI_Comm            comm;
601a914927fSHong Zhang   PetscMPIInt         size,rank,tagi,tagj,*len_si,*len_s,*len_ri,icompleted=0;
60242c57489SHong Zhang   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
60342c57489SHong Zhang   PetscInt            len,proc,*dnz,*onz,*owners;
60424ecddacSHong Zhang   PetscInt            nzi,*pti,*ptj;
60542c57489SHong Zhang   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
606ded4f561SHong Zhang   MPI_Request         *swaits,*rwaits;
607ded4f561SHong Zhang   MPI_Status          *sstatus,rstatus;
60842c57489SHong Zhang   Mat_Merge_SeqsToMPI *merge;
60977584fe6SHong Zhang   PetscInt            *api,*apj,*Jptr,apnz,*prmap=p->garray,pon,nspacedouble=0,j,ap_rmax=0;
610d16ca5f0SHong Zhang   PetscReal           afill=1.0,afill_tmp;
611a914927fSHong Zhang   PetscInt            rmax;
61242c57489SHong Zhang 
61342c57489SHong Zhang   PetscFunctionBegin;
614ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
61524ecddacSHong Zhang 
616c5af039cSHong Zhang   /* check if matrix local sizes are compatible */
617c5af039cSHong Zhang   if (A->rmap->rstart != P->rmap->rstart || A->rmap->rend != P->rmap->rend) {
618c5af039cSHong 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);
619c5af039cSHong Zhang   }
620c5af039cSHong Zhang   if (A->cmap->rstart != P->rmap->rstart || A->cmap->rend != P->rmap->rend) {
621c5af039cSHong 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);
622c5af039cSHong Zhang   }
623c5af039cSHong Zhang 
62483445d95SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
62583445d95SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
62683445d95SHong Zhang 
627f8487c73SHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
628b00a9115SJed Brown   ierr        = PetscNew(&ptap);CHKERRQ(ierr);
629b00a9115SJed Brown   ierr        = PetscNew(&merge);CHKERRQ(ierr);
6309f4155fbSHong Zhang   ptap->merge = merge;
631f8487c73SHong Zhang   ptap->reuse = MAT_INITIAL_MATRIX;
6326c6e00e0SHong Zhang 
6336c6e00e0SHong Zhang   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
634b7f45c76SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
635fe615159SHong Zhang 
6366c6e00e0SHong Zhang   /* get P_loc by taking all local rows of P */
637a1a4e84aSHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
6386c6e00e0SHong Zhang 
639a1a4e84aSHong Zhang   p_loc  = (Mat_SeqAIJ*)(ptap->P_loc)->data;
640a1a4e84aSHong Zhang   p_oth  = (Mat_SeqAIJ*)(ptap->P_oth)->data;
6416c6e00e0SHong Zhang   pi_loc = p_loc->i; pj_loc = p_loc->j;
6426c6e00e0SHong Zhang   pi_oth = p_oth->i; pj_oth = p_oth->j;
64342c57489SHong Zhang 
6442addb229SHong Zhang   /* (1) compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth (api,apj) */
6452addb229SHong Zhang   /*--------------------------------------------------------------------------*/
646854ce69bSBarry Smith   ierr   = PetscMalloc1(am+1,&api);CHKERRQ(ierr);
64782412ba7SHong Zhang   api[0] = 0;
64883445d95SHong Zhang 
64983445d95SHong Zhang   /* create and initialize a linked list */
650a914927fSHong Zhang   ierr = PetscLLCondensedCreate(pN,pN,&lnk,&lnkbt);CHKERRQ(ierr);
65183445d95SHong Zhang 
652a914927fSHong Zhang   /* Initial FreeSpace size is fill*(nnz(A) + nnz(P)) -OOM for ex56, np=8k on Intrepid! */
653d16ca5f0SHong Zhang   ierr = PetscFreeSpaceGet((PetscInt)(fill*(adi[am]+aoi[am]+pi_loc[pm])),&free_space);CHKERRQ(ierr);
6542205254eSKarl Rupp 
655f4a743e1SHong Zhang   current_space = free_space;
656f4a743e1SHong Zhang 
657f4a743e1SHong Zhang   for (i=0; i<am; i++) {
658f4a743e1SHong Zhang     /* diagonal portion of A */
659ded4f561SHong Zhang     nzi = adi[i+1] - adi[i];
66077584fe6SHong Zhang     aj  = ad->j + adi[i];
661ded4f561SHong Zhang     for (j=0; j<nzi; j++) {
66277584fe6SHong Zhang       row  = aj[j];
66382412ba7SHong Zhang       pnz  = pi_loc[row+1] - pi_loc[row];
664ded4f561SHong Zhang       Jptr = pj_loc + pi_loc[row];
66583445d95SHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
666a914927fSHong Zhang       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
667f4a743e1SHong Zhang     }
668f4a743e1SHong Zhang     /* off-diagonal portion of A */
669ded4f561SHong Zhang     nzi = aoi[i+1] - aoi[i];
67077584fe6SHong Zhang     aj  = ao->j + aoi[i];
671ded4f561SHong Zhang     for (j=0; j<nzi; j++) {
67277584fe6SHong Zhang       row  = aj[j];
67382412ba7SHong Zhang       pnz  = pi_oth[row+1] - pi_oth[row];
674ded4f561SHong Zhang       Jptr = pj_oth + pi_oth[row];
675a914927fSHong Zhang       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
676f4a743e1SHong Zhang     }
677a914927fSHong Zhang     apnz     = lnk[0];
67882412ba7SHong Zhang     api[i+1] = api[i] + apnz;
67977584fe6SHong Zhang     if (ap_rmax < apnz) ap_rmax = apnz;
680f4a743e1SHong Zhang 
68183445d95SHong Zhang     /* if free space is not available, double the total space in the list */
68282412ba7SHong Zhang     if (current_space->local_remaining<apnz) {
6832ba1acd5SHong Zhang       ierr = PetscFreeSpaceGet(apnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
684f2b054eeSHong Zhang       nspacedouble++;
685f4a743e1SHong Zhang     }
686f4a743e1SHong Zhang 
687f4a743e1SHong Zhang     /* Copy data into free space, then initialize lnk */
688a914927fSHong Zhang     ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
6892205254eSKarl Rupp 
69082412ba7SHong Zhang     current_space->array           += apnz;
69182412ba7SHong Zhang     current_space->local_used      += apnz;
69282412ba7SHong Zhang     current_space->local_remaining -= apnz;
693f4a743e1SHong Zhang   }
694a914927fSHong Zhang 
69582412ba7SHong Zhang   /* Allocate space for apj, initialize apj, and */
696f4a743e1SHong Zhang   /* destroy list of free space and other temporary array(s) */
697854ce69bSBarry Smith   ierr      = PetscMalloc1(api[am]+1,&apj);CHKERRQ(ierr);
69877584fe6SHong Zhang   ierr      = PetscFreeSpaceContiguous(&free_space,apj);CHKERRQ(ierr);
699118727c9SMark F. Adams   afill_tmp = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1);
700d16ca5f0SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
701f4a743e1SHong Zhang 
7022addb229SHong Zhang   /* (2) determine symbolic Co=(p->B)^T*AP - send to others (coi,coj)*/
7032addb229SHong Zhang   /*-----------------------------------------------------------------*/
704de0260b3SHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr);
70542c57489SHong Zhang 
706ded4f561SHong Zhang   /* then, compute symbolic Co = (p->B)^T*AP */
707d0f46423SBarry Smith   pon    = (p->B)->cmap->n; /* total num of rows to be sent to other processors
70883445d95SHong Zhang                          >= (num of nonzero rows of C_seq) - pn */
709854ce69bSBarry Smith   ierr   = PetscMalloc1(pon+1,&coi);CHKERRQ(ierr);
710de0260b3SHong Zhang   coi[0] = 0;
71142c57489SHong Zhang 
712d16ca5f0SHong Zhang   /* set initial free space to be fill*(nnz(p->B) + nnz(AP)) */
713d16ca5f0SHong Zhang   nnz           = fill*(poti[pon] + api[am]);
71422d28d08SBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
71542c57489SHong Zhang   current_space = free_space;
71642c57489SHong Zhang 
717de0260b3SHong Zhang   for (i=0; i<pon; i++) {
71882412ba7SHong Zhang     pnz = poti[i+1] - poti[i];
71977584fe6SHong Zhang     ptJ = potj + poti[i];
72077584fe6SHong Zhang     for (j=0; j<pnz; j++) {
72177584fe6SHong Zhang       row  = ptJ[j]; /* row of AP == col of Pot */
72282412ba7SHong Zhang       apnz = api[row+1] - api[row];
723ded4f561SHong Zhang       Jptr = apj + api[row];
72483445d95SHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
725a914927fSHong Zhang       ierr = PetscLLCondensedAddSorted(apnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
72642c57489SHong Zhang     }
727a914927fSHong Zhang     nnz = lnk[0];
72842c57489SHong Zhang 
72983445d95SHong Zhang     /* If free space is not available, double the total space in the list */
730ded4f561SHong Zhang     if (current_space->local_remaining<nnz) {
7314238b7adSHong Zhang       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
732d16ca5f0SHong Zhang       nspacedouble++;
73342c57489SHong Zhang     }
73442c57489SHong Zhang 
73542c57489SHong Zhang     /* Copy data into free space, and zero out denserows */
736a914927fSHong Zhang     ierr = PetscLLCondensedClean(pN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
7372205254eSKarl Rupp 
738ded4f561SHong Zhang     current_space->array           += nnz;
739ded4f561SHong Zhang     current_space->local_used      += nnz;
740ded4f561SHong Zhang     current_space->local_remaining -= nnz;
7412205254eSKarl Rupp 
742ded4f561SHong Zhang     coi[i+1] = coi[i] + nnz;
74342c57489SHong Zhang   }
7446b911d16SHong Zhang 
7456b911d16SHong Zhang   ierr      = PetscMalloc1(coi[pon],&coj);CHKERRQ(ierr);
746a1a86e44SBarry Smith   ierr      = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr);
747118727c9SMark F. Adams   afill_tmp = (PetscReal)coi[pon]/(poti[pon] + api[am]+1);
748d16ca5f0SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
749de0260b3SHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr);
75042c57489SHong Zhang 
7516b911d16SHong Zhang   /* (3) send j-array (coj) of Co to other processors */
7526b911d16SHong Zhang   /*--------------------------------------------------*/
7536b911d16SHong Zhang   ierr = PetscCalloc1(size,&merge->len_s);CHKERRQ(ierr);
7546b911d16SHong Zhang   len_s        = merge->len_s;
7556b911d16SHong Zhang   merge->nsend = 0;
7566b911d16SHong Zhang 
7576b911d16SHong Zhang 
75842c57489SHong Zhang   /* determine row ownership */
75926283091SBarry Smith   ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr);
7607a2fc3feSBarry Smith   merge->rowmap->n  = pn;
7617a2fc3feSBarry Smith   merge->rowmap->bs = 1;
7622205254eSKarl Rupp 
76326283091SBarry Smith   ierr   = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr);
7647a2fc3feSBarry Smith   owners = merge->rowmap->range;
76542c57489SHong Zhang 
76642c57489SHong Zhang   /* determine the number of messages to send, their lengths */
767dcca6d9dSJed Brown   ierr = PetscMalloc2(size,&len_si,size,&sstatus);CHKERRQ(ierr);
76883445d95SHong Zhang   ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
769854ce69bSBarry Smith   ierr = PetscMalloc1(size+2,&owners_co);CHKERRQ(ierr);
770de0260b3SHong Zhang 
77183445d95SHong Zhang   proc = 0;
772de0260b3SHong Zhang   for (i=0; i<pon; i++) {
773de0260b3SHong Zhang     while (prmap[i] >= owners[proc+1]) proc++;
7742addb229SHong Zhang     len_si[proc]++;               /* num of rows in Co(=Pt*AP) to be sent to [proc] */
775ee6e4b5bSHong Zhang     len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */
77683445d95SHong Zhang   }
777de0260b3SHong Zhang 
7782addb229SHong Zhang   len          = 0; /* max length of buf_si[], see (4) */
779de0260b3SHong Zhang   owners_co[0] = 0;
78042c57489SHong Zhang   for (proc=0; proc<size; proc++) {
781de0260b3SHong Zhang     owners_co[proc+1] = owners_co[proc] + len_si[proc];
782ee6e4b5bSHong Zhang     if (len_s[proc]) {
78342c57489SHong Zhang       merge->nsend++;
7842addb229SHong Zhang       len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */
78542c57489SHong Zhang       len         += len_si[proc];
78642c57489SHong Zhang     }
78742c57489SHong Zhang   }
78842c57489SHong Zhang 
789ded4f561SHong Zhang   /* determine the number and length of messages to receive for coi and coj  */
7900298fd71SBarry Smith   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&merge->nrecv);CHKERRQ(ierr);
79142c57489SHong Zhang   ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr);
79242c57489SHong Zhang 
793ded4f561SHong Zhang   /* post the Irecv and Isend of coj */
794529bc5dcSHong Zhang   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
795529bc5dcSHong Zhang   ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
796854ce69bSBarry Smith   ierr = PetscMalloc1(merge->nsend+1,&swaits);CHKERRQ(ierr);
79742c57489SHong Zhang   for (proc=0, k=0; proc<size; proc++) {
79842c57489SHong Zhang     if (!len_s[proc]) continue;
799de0260b3SHong Zhang     i    = owners_co[proc];
800529bc5dcSHong Zhang     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
80142c57489SHong Zhang     k++;
80242c57489SHong Zhang   }
80342c57489SHong Zhang 
804ded4f561SHong Zhang   /* receives and sends of coj are complete */
80577584fe6SHong Zhang   for (i=0; i<merge->nrecv; i++) {
806c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
807ded4f561SHong Zhang   }
808ded4f561SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
8090c468ba9SBarry Smith   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
81082412ba7SHong Zhang 
8116b911d16SHong Zhang   /* (4) send and recv coi */
8126b911d16SHong Zhang   /*-----------------------*/
813529bc5dcSHong Zhang   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
814529bc5dcSHong Zhang   ierr   = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
815854ce69bSBarry Smith   ierr   = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr);
81642c57489SHong Zhang   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
81742c57489SHong Zhang   for (proc=0,k=0; proc<size; proc++) {
81842c57489SHong Zhang     if (!len_s[proc]) continue;
81942c57489SHong Zhang     /* form outgoing message for i-structure:
82042c57489SHong Zhang          buf_si[0]:                 nrows to be sent
82142c57489SHong Zhang                [1:nrows]:           row index (global)
82242c57489SHong Zhang                [nrows+1:2*nrows+1]: i-structure index
82342c57489SHong Zhang     */
82442c57489SHong Zhang     /*-------------------------------------------*/
8252addb229SHong Zhang     nrows       = len_si[proc]/2 - 1; /* num of rows in Co to be sent to [proc] */
82642c57489SHong Zhang     buf_si_i    = buf_si + nrows+1;
82742c57489SHong Zhang     buf_si[0]   = nrows;
82842c57489SHong Zhang     buf_si_i[0] = 0;
82942c57489SHong Zhang     nrows       = 0;
830de0260b3SHong Zhang     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
831ded4f561SHong Zhang       nzi = coi[i+1] - coi[i];
832ded4f561SHong Zhang       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi;  /* i-structure */
833de0260b3SHong Zhang       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
83442c57489SHong Zhang       nrows++;
83542c57489SHong Zhang     }
836529bc5dcSHong Zhang     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
83742c57489SHong Zhang     k++;
83842c57489SHong Zhang     buf_si += len_si[proc];
83942c57489SHong Zhang   }
840ded4f561SHong Zhang   i = merge->nrecv;
841ded4f561SHong Zhang   while (i--) {
842c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
843ded4f561SHong Zhang   }
844ded4f561SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
8450c468ba9SBarry Smith   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
846a914927fSHong Zhang 
84724ecddacSHong Zhang   ierr = PetscFree2(len_si,sstatus);CHKERRQ(ierr);
84842c57489SHong Zhang   ierr = PetscFree(len_ri);CHKERRQ(ierr);
849ded4f561SHong Zhang   ierr = PetscFree(swaits);CHKERRQ(ierr);
85042c57489SHong Zhang   ierr = PetscFree(buf_s);CHKERRQ(ierr);
85142c57489SHong Zhang 
8526b911d16SHong Zhang   /* (5) compute the local portion of C (mpi mat) */
8536b911d16SHong Zhang   /*----------------------------------------------*/
854ded4f561SHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr);
855ded4f561SHong Zhang 
85624ecddacSHong Zhang   /* allocate pti array and free space for accumulating nonzero column info */
857854ce69bSBarry Smith   ierr   = PetscMalloc1(pn+1,&pti);CHKERRQ(ierr);
85824ecddacSHong Zhang   pti[0] = 0;
85942c57489SHong Zhang 
860d16ca5f0SHong Zhang   /* set initial free space to be fill*(nnz(P) + nnz(AP)) */
861d16ca5f0SHong Zhang   nnz           = fill*(pi_loc[pm] + api[am]);
86222d28d08SBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
86342c57489SHong Zhang   current_space = free_space;
86442c57489SHong Zhang 
865dcca6d9dSJed Brown   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
86642c57489SHong Zhang   for (k=0; k<merge->nrecv; k++) {
86742c57489SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
86842c57489SHong Zhang     nrows       = *buf_ri_k[k];
86942c57489SHong Zhang     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
87042c57489SHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
87142c57489SHong Zhang   }
87242c57489SHong Zhang   ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr);
87308cb4532SHong Zhang   rmax = 0;
87442c57489SHong Zhang   for (i=0; i<pn; i++) {
875ded4f561SHong Zhang     /* add pdt[i,:]*AP into lnk */
876ded4f561SHong Zhang     pnz = pdti[i+1] - pdti[i];
87777584fe6SHong Zhang     ptJ = pdtj + pdti[i];
87877584fe6SHong Zhang     for (j=0; j<pnz; j++) {
87977584fe6SHong Zhang       row  = ptJ[j];  /* row of AP == col of Pt */
880ded4f561SHong Zhang       apnz = api[row+1] - api[row];
881ded4f561SHong Zhang       Jptr = apj + api[row];
882ded4f561SHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
883a914927fSHong Zhang       ierr = PetscLLCondensedAddSorted(apnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
884ded4f561SHong Zhang     }
885a914927fSHong Zhang 
88642c57489SHong Zhang     /* add received col data into lnk */
88742c57489SHong Zhang     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
88842c57489SHong Zhang       if (i == *nextrow[k]) { /* i-th row */
889ded4f561SHong Zhang         nzi  = *(nextci[k]+1) - *nextci[k];
890ded4f561SHong Zhang         Jptr = buf_rj[k] + *nextci[k];
891a914927fSHong Zhang         ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr);
89242c57489SHong Zhang         nextrow[k]++; nextci[k]++;
89342c57489SHong Zhang       }
89442c57489SHong Zhang     }
895a914927fSHong Zhang     nnz = lnk[0];
89642c57489SHong Zhang 
89742c57489SHong Zhang     /* if free space is not available, make more free space */
898ded4f561SHong Zhang     if (current_space->local_remaining<nnz) {
8994238b7adSHong Zhang       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
900d16ca5f0SHong Zhang       nspacedouble++;
90142c57489SHong Zhang     }
90242c57489SHong Zhang     /* copy data into free space, then initialize lnk */
903a914927fSHong Zhang     ierr = PetscLLCondensedClean(pN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
904ded4f561SHong Zhang     ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr);
9052205254eSKarl Rupp 
906ded4f561SHong Zhang     current_space->array           += nnz;
907ded4f561SHong Zhang     current_space->local_used      += nnz;
908ded4f561SHong Zhang     current_space->local_remaining -= nnz;
9092205254eSKarl Rupp 
91024ecddacSHong Zhang     pti[i+1] = pti[i] + nnz;
91108cb4532SHong Zhang     if (nnz > rmax) rmax = nnz;
91242c57489SHong Zhang   }
913ded4f561SHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr);
9140572522cSBarry Smith   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
91542c57489SHong Zhang 
916854ce69bSBarry Smith   ierr      = PetscMalloc1(pti[pn]+1,&ptj);CHKERRQ(ierr);
91724ecddacSHong Zhang   ierr      = PetscFreeSpaceContiguous(&free_space,ptj);CHKERRQ(ierr);
91824ecddacSHong Zhang   afill_tmp = (PetscReal)pti[pn]/(pi_loc[pm] + api[am]+1);
919d16ca5f0SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
92042c57489SHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
92142c57489SHong Zhang 
9226b911d16SHong Zhang   /* (6) create symbolic parallel matrix Cmpi */
9236b911d16SHong Zhang   /*------------------------------------------*/
92477584fe6SHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
92577584fe6SHong Zhang   ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
92633d57670SJed Brown   ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr);
92777584fe6SHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
92877584fe6SHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
92942c57489SHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
930a2f3521dSMark F. Adams 
931b091e509SHong Zhang   merge->bi        = pti;      /* Cseq->i */
932b091e509SHong Zhang   merge->bj        = ptj;      /* Cseq->j */
933b091e509SHong Zhang   merge->coi       = coi;      /* Co->i   */
934b091e509SHong Zhang   merge->coj       = coj;      /* Co->j   */
93542c57489SHong Zhang   merge->buf_ri    = buf_ri;
93642c57489SHong Zhang   merge->buf_rj    = buf_rj;
937de0260b3SHong Zhang   merge->owners_co = owners_co;
93877584fe6SHong Zhang   merge->destroy   = Cmpi->ops->destroy;
93977584fe6SHong Zhang   merge->duplicate = Cmpi->ops->duplicate;
940cc6cb787SHong Zhang 
94177584fe6SHong Zhang   /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */
942a914927fSHong Zhang   Cmpi->assembled      = PETSC_FALSE;
94377584fe6SHong Zhang   Cmpi->ops->destroy   = MatDestroy_MPIAIJ_PtAP;
94477584fe6SHong Zhang   Cmpi->ops->duplicate = MatDuplicate_MPIAIJ_MatPtAP;
94542c57489SHong Zhang 
94677584fe6SHong Zhang   /* attach the supporting struct to Cmpi for reuse */
94777584fe6SHong Zhang   c           = (Mat_MPIAIJ*)Cmpi->data;
948f8487c73SHong Zhang   c->ptap     = ptap;
94977584fe6SHong Zhang   ptap->api   = api;
95077584fe6SHong Zhang   ptap->apj   = apj;
951d6ab1dc0SHong Zhang   ptap->rmax  = ap_rmax;
95277584fe6SHong Zhang   *C          = Cmpi;
953414894bdSHong Zhang 
954414894bdSHong Zhang   /* flag 'scalable' determines which implementations to be used:
955414894bdSHong Zhang        0: do dense axpy in MatPtAPNumeric() - fast, but requires storage of a nonscalable dense array apa;
956414894bdSHong Zhang        1: do sparse axpy in MatPtAPNumeric() - might slow, uses a sparse array apa */
957414894bdSHong Zhang   /* set default scalable */
9589516a85cSHong Zhang   ptap->scalable = PETSC_FALSE; //PETSC_TRUE;
9592205254eSKarl Rupp 
9600298fd71SBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)Cmpi)->prefix,"-matptap_scalable",&ptap->scalable,NULL);CHKERRQ(ierr);
961414894bdSHong Zhang   if (!ptap->scalable) {  /* Do dense axpy */
9621795a4d1SJed Brown     ierr = PetscCalloc1(pN,&ptap->apa);CHKERRQ(ierr);
963414894bdSHong Zhang   } else {
9641795a4d1SJed Brown     ierr = PetscCalloc1(ap_rmax+1,&ptap->apa);CHKERRQ(ierr);
965414894bdSHong Zhang   }
966414894bdSHong Zhang 
967f2b054eeSHong Zhang #if defined(PETSC_USE_INFO)
96824ecddacSHong Zhang   if (pti[pn] != 0) {
96957622a8eSBarry Smith     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
97057622a8eSBarry Smith     ierr = PetscInfo1(Cmpi,"Use MatPtAP(A,P,MatReuse,%g,&C) for best performance.\n",(double)afill);CHKERRQ(ierr);
971f2b054eeSHong Zhang   } else {
97277584fe6SHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
973f2b054eeSHong Zhang   }
974f2b054eeSHong Zhang #endif
97542c57489SHong Zhang   PetscFunctionReturn(0);
97642c57489SHong Zhang }
97742c57489SHong Zhang 
97842c57489SHong Zhang #undef __FUNCT__
97942c57489SHong Zhang #define __FUNCT__ "MatPtAPNumeric_MPIAIJ_MPIAIJ"
98042c57489SHong Zhang PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C)
98142c57489SHong Zhang {
98242c57489SHong Zhang   PetscErrorCode      ierr;
983f8487c73SHong Zhang   Mat_MPIAIJ          *a =(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
98442c57489SHong Zhang   Mat_SeqAIJ          *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
985de0260b3SHong Zhang   Mat_SeqAIJ          *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data;
98642c57489SHong Zhang   Mat_SeqAIJ          *p_loc,*p_oth;
987f8487c73SHong Zhang   Mat_PtAPMPI         *ptap;
9889f4155fbSHong Zhang   Mat_Merge_SeqsToMPI *merge;
9891d633602SHong Zhang   PetscInt            *adi=ad->i,*aoi=ao->i,*adj,*aoj,*apJ,nextp;
99082412ba7SHong Zhang   PetscInt            *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pJ,*pj;
99182412ba7SHong Zhang   PetscInt            i,j,k,anz,pnz,apnz,nextap,row,*cj;
992e900a757SHong Zhang   MatScalar           *ada,*aoa,*apa,*pa,*ca,*pa_loc,*pa_oth,valtmp;
993d0f46423SBarry Smith   PetscInt            am  =A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n;
994ce94432eSBarry Smith   MPI_Comm            comm;
99542c57489SHong Zhang   PetscMPIInt         size,rank,taga,*len_s;
99682412ba7SHong Zhang   PetscInt            *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci;
99742c57489SHong Zhang   PetscInt            **buf_ri,**buf_rj;
99850f5bf86SHong Zhang   PetscInt            cnz=0,*bj_i,*bi,*bj,bnz,nextcj;  /* bi,bj,ba: local array of C(mpi mat) */
99942c57489SHong Zhang   MPI_Request         *s_waits,*r_waits;
100042c57489SHong Zhang   MPI_Status          *status;
100182412ba7SHong Zhang   MatScalar           **abuf_r,*ba_i,*pA,*coa,*ba;
100282412ba7SHong Zhang   PetscInt            *api,*apj,*coi,*coj;
1003d0f46423SBarry Smith   PetscInt            *poJ=po->j,*pdJ=pd->j,pcstart=P->cmap->rstart,pcend=P->cmap->rend;
10046ce94e8fSHong Zhang   PetscBool           scalable;
100538571addSHong Zhang #if defined(PTAP_PROFILE)
1006e497d3c8SHong Zhang   PetscLogDouble t0,t1,t2,eP,t3,t4,et2_AP=0.0,ePtAP=0.0,t2_0,t2_1,t2_2;
100738571addSHong Zhang #endif
100842c57489SHong Zhang 
100942c57489SHong Zhang   PetscFunctionBegin;
1010ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr);
101138571addSHong Zhang #if defined(PTAP_PROFILE)
10128563dfccSBarry Smith   ierr = PetscTime(&t0);CHKERRQ(ierr);
101338571addSHong Zhang #endif
101442c57489SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
101542c57489SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
101642c57489SHong Zhang 
1017f8487c73SHong Zhang   ptap = c->ptap;
1018ce94432eSBarry 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");
1019f8487c73SHong Zhang   merge    = ptap->merge;
1020414894bdSHong Zhang   apa      = ptap->apa;
10216ce94e8fSHong Zhang   scalable = ptap->scalable;
10226c6e00e0SHong Zhang 
1023a1a4e84aSHong Zhang   /* 1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
10246b911d16SHong Zhang   /*-----------------------------------------------------*/
1025f8487c73SHong Zhang   if (ptap->reuse == MAT_INITIAL_MATRIX) {
10269f4155fbSHong Zhang     /* P_oth and P_loc are obtained in MatPtASymbolic(), skip calling MatGetBrowsOfAoCols() and MatMPIAIJGetLocalMat() */
1027f8487c73SHong Zhang     ptap->reuse = MAT_REUSE_MATRIX;
10286c6e00e0SHong Zhang   } else { /* update numerical values of P_oth and P_loc */
1029b7f45c76SHong Zhang     ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
1030a1a4e84aSHong Zhang     ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
10316c6e00e0SHong Zhang   }
103238571addSHong Zhang #if defined(PTAP_PROFILE)
10338563dfccSBarry Smith   ierr = PetscTime(&t1);CHKERRQ(ierr);
103405d62848SHong Zhang   eP = t1-t0;
103538571addSHong Zhang #endif
103605d62848SHong Zhang   /*
103705d62848SHong Zhang   printf("[%d] Ad: %d, %d; Ao: %d, %d; P_loc: %d, %d; P_oth %d, %d;\n",rank,
103805d62848SHong Zhang          a->A->rmap->N,a->A->cmap->N,a->B->rmap->N,a->B->cmap->N,
103905d62848SHong Zhang          ptap->P_loc->rmap->N,ptap->P_loc->cmap->N,
104005d62848SHong Zhang          ptap->P_oth->rmap->N,ptap->P_oth->cmap->N);
104105d62848SHong Zhang    */
1042f8487c73SHong Zhang 
1043f9473cf7SHong Zhang   /* 2) compute numeric C_seq = P_loc^T*A_loc*P - dominating part */
1044f9473cf7SHong Zhang   /*--------------------------------------------------------------*/
104542c57489SHong Zhang   /* get data from symbolic products */
1046a1a4e84aSHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
1047a1a4e84aSHong Zhang   p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
1048a9d06231SHong Zhang   pi_loc=p_loc->i; pj_loc=p_loc->j; pJ=pj_loc; pa_loc=p_loc->a;
104942c57489SHong Zhang   pi_oth=p_oth->i; pj_oth=p_oth->j; pa_oth=p_oth->a;
105042c57489SHong Zhang 
1051de0260b3SHong Zhang   coi  = merge->coi; coj = merge->coj;
10521795a4d1SJed Brown   ierr = PetscCalloc1(coi[pon]+1,&coa);CHKERRQ(ierr);
1053de0260b3SHong Zhang 
1054de0260b3SHong Zhang   bi     = merge->bi; bj = merge->bj;
10557a2fc3feSBarry Smith   owners = merge->rowmap->range;
10561795a4d1SJed Brown   ierr   = PetscCalloc1(bi[cm]+1,&ba);CHKERRQ(ierr);  /* ba: Cseq->a */
105742c57489SHong Zhang 
1058a1a4e84aSHong Zhang   api = ptap->api; apj = ptap->apj;
1059d9824c15SHong Zhang 
10609516a85cSHong Zhang   if (!scalable) { /* Do dense axpy on apa (length of pN, stores A[i,:]*P) - nonscalable, but faster (could take 1/3 scalable time) */
1061b38d1a2dSJed Brown     ierr = PetscInfo(C,"Using non-scalable dense axpy\n");CHKERRQ(ierr);
106205d62848SHong Zhang #if 0
10630d3441aeSHong Zhang     /* ------ 10x slower -------------- */
10640d3441aeSHong Zhang     /*==================================*/
10659516a85cSHong Zhang     Mat         R = ptap->R;
10669516a85cSHong Zhang     Mat_SeqAIJ  *r = (Mat_SeqAIJ*)R->data;
10679516a85cSHong Zhang     PetscInt    *ri=r->i,*rj=r->j,rnz,arow,l,prow,pcol,pN=P->cmap->N;
10689516a85cSHong Zhang     PetscScalar *ra=r->a,tmp,cdense[pN];
10699516a85cSHong Zhang 
10709516a85cSHong Zhang     ierr = PetscMemzero(cdense,pN*sizeof(PetscScalar));CHKERRQ(ierr);
10719516a85cSHong Zhang     for (i=0; i<cm; i++) { /* each row of C or R */
10729516a85cSHong Zhang       rnz = ri[i+1] - ri[i];
10739516a85cSHong Zhang 
10749516a85cSHong Zhang       for (j=0; j<rnz; j++) { /* each nz of R */
10759516a85cSHong Zhang         arow = rj[ri[i] + j];
10769516a85cSHong Zhang 
10779516a85cSHong Zhang         /* diagonal portion of A */
10789516a85cSHong Zhang         anz  = ad->i[arow+1] - ad->i[arow];
10799516a85cSHong Zhang         for (k=0; k<anz; k++) { /* each nz of Ad */
10809516a85cSHong Zhang           tmp  = ra[ri[i] + j]*ad->a[ad->i[arow] + k];
10819516a85cSHong Zhang           prow = ad->j[ad->i[arow] + k];
10829516a85cSHong Zhang           pnz  = pi_loc[prow+1] - pi_loc[prow];
10839516a85cSHong Zhang 
10849516a85cSHong Zhang           for (l=0; l<pnz; l++) { /* each nz of P_loc */
10859516a85cSHong Zhang             pcol = pj_loc[pi_loc[prow] + l];
10869516a85cSHong Zhang             cdense[pcol] += tmp*pa_loc[pi_loc[prow] + l];
10879516a85cSHong Zhang           }
10889516a85cSHong Zhang         }
10899516a85cSHong Zhang 
10909516a85cSHong Zhang         /* off-diagonal portion of A */
10919516a85cSHong Zhang         anz  = ao->i[arow+1] - ao->i[arow];
10929516a85cSHong Zhang         for (k=0; k<anz; k++) { /* each nz of Ao */
10939516a85cSHong Zhang           tmp  = ra[ri[i] + j]*ao->a[ao->i[arow] + k];
10949516a85cSHong Zhang           prow = ao->j[ao->i[arow] + k];
10959516a85cSHong Zhang           pnz  = pi_oth[prow+1] - pi_oth[prow];
10969516a85cSHong Zhang 
10979516a85cSHong Zhang           for (l=0; l<pnz; l++) { /* each nz of P_oth */
10989516a85cSHong Zhang             pcol = pj_oth[pi_oth[prow] + l];
10999516a85cSHong Zhang             cdense[pcol] += tmp*pa_oth[pi_oth[prow] + l];
11009516a85cSHong Zhang           }
11019516a85cSHong Zhang         }
11029516a85cSHong Zhang 
11039516a85cSHong Zhang       } //for (j=0; j<rnz; j++)
11049516a85cSHong Zhang 
11059516a85cSHong Zhang       /* copy cdense[] into ca; zero cdense[] */
11069516a85cSHong Zhang       cnz = bi[i+1] - bi[i];
11079516a85cSHong Zhang       cj  = bj + bi[i];
11089516a85cSHong Zhang       ca  = ba + bi[i];
11099516a85cSHong Zhang       for (j=0; j<cnz; j++) {
11109516a85cSHong Zhang         ca[j] += cdense[cj[j]];
11119516a85cSHong Zhang         cdense[cj[j]] = 0.0;
11129516a85cSHong Zhang       }
11139516a85cSHong Zhang #if 0
11149516a85cSHong Zhang       if (rank == 0) {
11159516a85cSHong Zhang         printf("[%d] row %d: ",rank,i);
11169516a85cSHong Zhang         for (j=0; j<pN; j++) printf(" %g,",cdense[j]);
11179516a85cSHong Zhang         printf("\n");
11189516a85cSHong Zhang       }
11199516a85cSHong Zhang       for (j=0; j<pN; j++) cdense[j]=0.0; // zero cdnese[]
11209516a85cSHong Zhang #endif
11219516a85cSHong Zhang     } //for (i=0; i<cm; i++) {
112205d62848SHong Zhang #endif
11239516a85cSHong Zhang 
11249516a85cSHong Zhang     //==========================================
11259516a85cSHong Zhang 
112605d62848SHong Zhang     ierr = PetscTime(&t1);CHKERRQ(ierr);
1127a9d06231SHong Zhang     for (i=0; i<am; i++) {
1128b091e509SHong Zhang #if defined(PTAP_PROFILE)
11298563dfccSBarry Smith       ierr = PetscTime(&t2_0);CHKERRQ(ierr);
1130b091e509SHong Zhang #endif
1131a9d06231SHong Zhang       /* 2-a) form i-th sparse row of A_loc*P = Ad*P_loc + Ao*P_oth */
1132a9d06231SHong Zhang       /*------------------------------------------------------------*/
1133a9d06231SHong Zhang       apJ = apj + api[i];
1134a9d06231SHong Zhang 
1135a9d06231SHong Zhang       /* diagonal portion of A */
1136a9d06231SHong Zhang       anz = adi[i+1] - adi[i];
1137a9d06231SHong Zhang       adj = ad->j + adi[i];
1138a9d06231SHong Zhang       ada = ad->a + adi[i];
1139a9d06231SHong Zhang       for (j=0; j<anz; j++) {
1140a9d06231SHong Zhang         row = adj[j];
1141a9d06231SHong Zhang         pnz = pi_loc[row+1] - pi_loc[row];
1142a9d06231SHong Zhang         pj  = pj_loc + pi_loc[row];
1143a9d06231SHong Zhang         pa  = pa_loc + pi_loc[row];
1144a9d06231SHong Zhang 
1145a9d06231SHong Zhang         /* perform dense axpy */
1146e900a757SHong Zhang         valtmp = ada[j];
1147a9d06231SHong Zhang         for (k=0; k<pnz; k++) {
1148e900a757SHong Zhang           apa[pj[k]] += valtmp*pa[k];
1149a9d06231SHong Zhang         }
1150a9d06231SHong Zhang         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
1151a9d06231SHong Zhang       }
1152a9d06231SHong Zhang 
1153a9d06231SHong Zhang       /* off-diagonal portion of A */
1154a9d06231SHong Zhang       anz = aoi[i+1] - aoi[i];
1155a9d06231SHong Zhang       aoj = ao->j + aoi[i];
1156a9d06231SHong Zhang       aoa = ao->a + aoi[i];
1157a9d06231SHong Zhang       for (j=0; j<anz; j++) {
1158a9d06231SHong Zhang         row = aoj[j];
1159a9d06231SHong Zhang         pnz = pi_oth[row+1] - pi_oth[row];
1160a9d06231SHong Zhang         pj  = pj_oth + pi_oth[row];
1161a9d06231SHong Zhang         pa  = pa_oth + pi_oth[row];
1162a9d06231SHong Zhang 
1163a9d06231SHong Zhang         /* perform dense axpy */
1164e900a757SHong Zhang         valtmp = aoa[j];
1165a9d06231SHong Zhang         for (k=0; k<pnz; k++) {
1166e900a757SHong Zhang           apa[pj[k]] += valtmp*pa[k];
1167a9d06231SHong Zhang         }
1168a9d06231SHong Zhang         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
1169a9d06231SHong Zhang       }
1170b091e509SHong Zhang #if defined(PTAP_PROFILE)
11718563dfccSBarry Smith       ierr    = PetscTime(&t2_1);CHKERRQ(ierr);
1172b091e509SHong Zhang       et2_AP += t2_1 - t2_0;
1173b091e509SHong Zhang #endif
1174a9d06231SHong Zhang 
1175a9d06231SHong Zhang       /* 2-b) Compute Cseq = P_loc[i,:]^T*AP[i,:] using outer product */
1176a9d06231SHong Zhang       /*--------------------------------------------------------------*/
1177a9d06231SHong Zhang       apnz = api[i+1] - api[i];
1178a9d06231SHong Zhang       /* put the value into Co=(p->B)^T*AP (off-diagonal part, send to others) */
1179a9d06231SHong Zhang       pnz = po->i[i+1] - po->i[i];
1180e900a757SHong Zhang       poJ = po->j + po->i[i];
1181a9d06231SHong Zhang       pA  = po->a + po->i[i];
1182a9d06231SHong Zhang       for (j=0; j<pnz; j++) {
1183e900a757SHong Zhang         row = poJ[j];
1184e900a757SHong Zhang         cnz = coi[row+1] - coi[row];
1185e900a757SHong Zhang         cj  = coj + coi[row];
1186e900a757SHong Zhang         ca  = coa + coi[row];
1187a9d06231SHong Zhang         /* perform dense axpy */
1188e900a757SHong Zhang         valtmp = pA[j];
1189a9d06231SHong Zhang         for (k=0; k<cnz; k++) {
1190e900a757SHong Zhang           ca[k] += valtmp*apa[cj[k]];
1191a9d06231SHong Zhang         }
1192a9d06231SHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1193a9d06231SHong Zhang       }
119405d62848SHong Zhang #if 1
1195a9d06231SHong Zhang       /* put the value into Cd (diagonal part) */
1196a9d06231SHong Zhang       pnz = pd->i[i+1] - pd->i[i];
1197e900a757SHong Zhang       pdJ = pd->j + pd->i[i];
1198a9d06231SHong Zhang       pA  = pd->a + pd->i[i];
1199a9d06231SHong Zhang       for (j=0; j<pnz; j++) {
1200e900a757SHong Zhang         row = pdJ[j];
1201e900a757SHong Zhang         cnz = bi[row+1] - bi[row];
1202e900a757SHong Zhang         cj  = bj + bi[row];
1203e900a757SHong Zhang         ca  = ba + bi[row];
1204a9d06231SHong Zhang         /* perform dense axpy */
1205e900a757SHong Zhang         valtmp = pA[j];
1206a9d06231SHong Zhang         for (k=0; k<cnz; k++) {
1207e900a757SHong Zhang           ca[k] += valtmp*apa[cj[k]];
1208a9d06231SHong Zhang         }
1209a9d06231SHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1210a9d06231SHong Zhang       }
12119516a85cSHong Zhang #endif
1212a9d06231SHong Zhang       /* zero the current row of A*P */
1213a9d06231SHong Zhang       for (k=0; k<apnz; k++) apa[apJ[k]] = 0.0;
1214b091e509SHong Zhang #if defined(PTAP_PROFILE)
12158563dfccSBarry Smith       ierr      = PetscTime(&t2_2);CHKERRQ(ierr);
121605d62848SHong Zhang       ePtAP += t2_2 - t2_1;
1217b091e509SHong Zhang #endif
1218a9d06231SHong Zhang     }
12199516a85cSHong Zhang 
12209516a85cSHong Zhang     if (rank == 100) {
12219516a85cSHong Zhang     for (row=0; row<cm; row++) {
12229516a85cSHong Zhang       printf("[%d] row %d: ",rank,row);
12239516a85cSHong Zhang       cnz = bi[row+1] - bi[row];
12249516a85cSHong Zhang       for (j=0; j<cnz; j++) printf(" %g,",ba[bi[row]+j]);
12259516a85cSHong Zhang       printf("\n");
12269516a85cSHong Zhang     }
12279516a85cSHong Zhang     }
12289516a85cSHong Zhang 
122938571addSHong Zhang   } else { /* Do sparse axpy on apa (length of ap_rmax, stores A[i,:]*P) - scalable, but slower */
1230b38d1a2dSJed Brown     ierr = PetscInfo(C,"Using scalable sparse axpy\n");CHKERRQ(ierr);
123138571addSHong Zhang     /*-----------------------------------------------------------------------------------------*/
1232a9d06231SHong Zhang     pA=pa_loc;
123342c57489SHong Zhang     for (i=0; i<am; i++) {
1234b091e509SHong Zhang #if defined(PTAP_PROFILE)
12358563dfccSBarry Smith       ierr = PetscTime(&t2_0);CHKERRQ(ierr);
1236b091e509SHong Zhang #endif
1237f4a743e1SHong Zhang       /* form i-th sparse row of A*P */
123882412ba7SHong Zhang       apnz = api[i+1] - api[i];
123982412ba7SHong Zhang       apJ  = apj + api[i];
124042c57489SHong Zhang       /* diagonal portion of A */
124182412ba7SHong Zhang       anz = adi[i+1] - adi[i];
12421d633602SHong Zhang       adj = ad->j + adi[i];
12431d633602SHong Zhang       ada = ad->a + adi[i];
124482412ba7SHong Zhang       for (j=0; j<anz; j++) {
12451d633602SHong Zhang         row    = adj[j];
124682412ba7SHong Zhang         pnz    = pi_loc[row+1] - pi_loc[row];
124782412ba7SHong Zhang         pj     = pj_loc + pi_loc[row];
124882412ba7SHong Zhang         pa     = pa_loc + pi_loc[row];
1249e900a757SHong Zhang         valtmp = ada[j];
1250f4a743e1SHong Zhang         nextp  = 0;
125182412ba7SHong Zhang         for (k=0; nextp<pnz; k++) {
125282412ba7SHong Zhang           if (apJ[k] == pj[nextp]) { /* col of AP == col of P */
1253e900a757SHong Zhang             apa[k] += valtmp*pa[nextp++];
125442c57489SHong Zhang           }
125542c57489SHong Zhang         }
1256dc0b31edSSatish Balay         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
125742c57489SHong Zhang       }
125842c57489SHong Zhang       /* off-diagonal portion of A */
125982412ba7SHong Zhang       anz = aoi[i+1] - aoi[i];
12601d633602SHong Zhang       aoj = ao->j + aoi[i];
12611d633602SHong Zhang       aoa = ao->a + aoi[i];
126282412ba7SHong Zhang       for (j=0; j<anz; j++) {
12631d633602SHong Zhang         row    = aoj[j];
126482412ba7SHong Zhang         pnz    = pi_oth[row+1] - pi_oth[row];
126582412ba7SHong Zhang         pj     = pj_oth + pi_oth[row];
126682412ba7SHong Zhang         pa     = pa_oth + pi_oth[row];
1267e900a757SHong Zhang         valtmp = aoa[j];
1268f4a743e1SHong Zhang         nextp  = 0;
126982412ba7SHong Zhang         for (k=0; nextp<pnz; k++) {
127082412ba7SHong Zhang           if (apJ[k] == pj[nextp]) { /* col of AP == col of P */
1271e900a757SHong Zhang             apa[k] += valtmp*pa[nextp++];
127242c57489SHong Zhang           }
127342c57489SHong Zhang         }
1274dc0b31edSSatish Balay         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
127542c57489SHong Zhang       }
1276b091e509SHong Zhang #if defined(PTAP_PROFILE)
12778563dfccSBarry Smith       ierr    = PetscTime(&t2_1);CHKERRQ(ierr);
1278b091e509SHong Zhang       et2_AP += t2_1 - t2_0;
1279b091e509SHong Zhang #endif
128042c57489SHong Zhang 
1281a9d06231SHong Zhang       /* 2-b) Compute Cseq = P_loc[i,:]^T*AP[i,:] using outer product */
1282a9d06231SHong Zhang       /*--------------------------------------------------------------*/
128382412ba7SHong Zhang       pnz = pi_loc[i+1] - pi_loc[i];
1284f9473cf7SHong Zhang       pJ  = pj_loc + pi_loc[i];
128582412ba7SHong Zhang       for (j=0; j<pnz; j++) {
128642c57489SHong Zhang         nextap = 0;
1287f9473cf7SHong Zhang         row    = pJ[j]; /* global index */
128882412ba7SHong Zhang         if (row < pcstart || row >=pcend) { /* put the value into Co */
1289e900a757SHong Zhang           row = *poJ;
1290e900a757SHong Zhang           cj  = coj + coi[row];
1291e900a757SHong Zhang           ca  = coa + coi[row]; poJ++;
129282412ba7SHong Zhang         } else {                            /* put the value into Cd */
1293e900a757SHong Zhang           row = *pdJ;
1294e900a757SHong Zhang           cj  = bj + bi[row];
1295e900a757SHong Zhang           ca  = ba + bi[row]; pdJ++;
129642c57489SHong Zhang         }
1297e900a757SHong Zhang         valtmp = pA[j];
129882412ba7SHong Zhang         for (k=0; nextap<apnz; k++) {
1299e900a757SHong Zhang           if (cj[k]==apJ[nextap]) ca[k] += valtmp*apa[nextap++];
130042c57489SHong Zhang         }
1301dc0b31edSSatish Balay         ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr);
1302de0260b3SHong Zhang       }
13030496f32cSHong Zhang       pA += pnz;
130442c57489SHong Zhang       /* zero the current row info for A*P */
130582412ba7SHong Zhang       ierr = PetscMemzero(apa,apnz*sizeof(MatScalar));CHKERRQ(ierr);
1306b091e509SHong Zhang #if defined(PTAP_PROFILE)
13078563dfccSBarry Smith       ierr      = PetscTime(&t2_2);CHKERRQ(ierr);
130805d62848SHong Zhang       ePtAP += t2_2 - t2_1;
1309b091e509SHong Zhang #endif
131042c57489SHong Zhang     }
1311d9824c15SHong Zhang   }
131238571addSHong Zhang #if defined(PTAP_PROFILE)
13138563dfccSBarry Smith   ierr = PetscTime(&t2);CHKERRQ(ierr);
131438571addSHong Zhang #endif
131542c57489SHong Zhang 
1316a9d06231SHong Zhang   /* 3) send and recv matrix values coa */
1317a9d06231SHong Zhang   /*------------------------------------*/
131842c57489SHong Zhang   buf_ri = merge->buf_ri;
131942c57489SHong Zhang   buf_rj = merge->buf_rj;
132042c57489SHong Zhang   len_s  = merge->len_s;
1321fc42d0c8SSatish Balay   ierr   = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr);
132242c57489SHong Zhang   ierr   = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr);
132342c57489SHong Zhang 
1324dcca6d9dSJed Brown   ierr = PetscMalloc2(merge->nsend+1,&s_waits,size,&status);CHKERRQ(ierr);
132542c57489SHong Zhang   for (proc=0,k=0; proc<size; proc++) {
132642c57489SHong Zhang     if (!len_s[proc]) continue;
1327de0260b3SHong Zhang     i    = merge->owners_co[proc];
1328de0260b3SHong Zhang     ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr);
132942c57489SHong Zhang     k++;
133042c57489SHong Zhang   }
13310c468ba9SBarry Smith   if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);}
13320c468ba9SBarry Smith   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);}
133342c57489SHong Zhang 
1334a9d06231SHong Zhang   ierr = PetscFree2(s_waits,status);CHKERRQ(ierr);
133542c57489SHong Zhang   ierr = PetscFree(r_waits);CHKERRQ(ierr);
133682412ba7SHong Zhang   ierr = PetscFree(coa);CHKERRQ(ierr);
133738571addSHong Zhang #if defined(PTAP_PROFILE)
13388563dfccSBarry Smith   ierr = PetscTime(&t3);CHKERRQ(ierr);
133938571addSHong Zhang #endif
134042c57489SHong Zhang 
1341a9d06231SHong Zhang   /* 4) insert local Cseq and received values into Cmpi */
1342a9d06231SHong Zhang   /*------------------------------------------------------*/
1343dcca6d9dSJed Brown   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
134442c57489SHong Zhang   for (k=0; k<merge->nrecv; k++) {
134542c57489SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
134642c57489SHong Zhang     nrows       = *(buf_ri_k[k]);
134742c57489SHong Zhang     nextrow[k]  = buf_ri_k[k]+1;  /* next row number of k-th recved i-structure */
134882412ba7SHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
134942c57489SHong Zhang   }
135042c57489SHong Zhang 
1351de0260b3SHong Zhang   for (i=0; i<cm; i++) {
135282412ba7SHong Zhang     row  = owners[rank] + i; /* global row index of C_seq */
135342c57489SHong Zhang     bj_i = bj + bi[i];  /* col indices of the i-th row of C */
1354de0260b3SHong Zhang     ba_i = ba + bi[i];
135582412ba7SHong Zhang     bnz  = bi[i+1] - bi[i];
135642c57489SHong Zhang     /* add received vals into ba */
135742c57489SHong Zhang     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
135842c57489SHong Zhang       /* i-th row */
135942c57489SHong Zhang       if (i == *nextrow[k]) {
136082412ba7SHong Zhang         cnz    = *(nextci[k]+1) - *nextci[k];
136182412ba7SHong Zhang         cj     = buf_rj[k] + *(nextci[k]);
136282412ba7SHong Zhang         ca     = abuf_r[k] + *(nextci[k]);
136382412ba7SHong Zhang         nextcj = 0;
136482412ba7SHong Zhang         for (j=0; nextcj<cnz; j++) {
136582412ba7SHong Zhang           if (bj_i[j] == cj[nextcj]) { /* bcol == ccol */
136682412ba7SHong Zhang             ba_i[j] += ca[nextcj++];
136742c57489SHong Zhang           }
136842c57489SHong Zhang         }
136982412ba7SHong Zhang         nextrow[k]++; nextci[k]++;
1370c5af039cSHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
137142c57489SHong Zhang       }
137242c57489SHong Zhang     }
137382412ba7SHong Zhang     ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr);
137442c57489SHong Zhang   }
137542c57489SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
137642c57489SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
137742c57489SHong Zhang 
1378de0260b3SHong Zhang   ierr = PetscFree(ba);CHKERRQ(ierr);
1379c05d87d6SBarry Smith   ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr);
138042c57489SHong Zhang   ierr = PetscFree(abuf_r);CHKERRQ(ierr);
13810572522cSBarry Smith   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
138238571addSHong Zhang #if defined(PTAP_PROFILE)
13838563dfccSBarry Smith   ierr = PetscTime(&t4);CHKERRQ(ierr);
13849516a85cSHong Zhang   if (rank==1) {
138505d62848SHong 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);
13869516a85cSHong Zhang   }
138738571addSHong Zhang #endif
138842c57489SHong Zhang   PetscFunctionReturn(0);
138942c57489SHong Zhang }
1390