xref: /petsc/src/mat/impls/aij/mpi/mpiptap.c (revision de817e96340c6676a88709a354a002e2c99a76b3)
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;
1180d3441aeSHong Zhang 
1190d3441aeSHong Zhang   PetscFunctionBegin;
1202259aa2eSHong Zhang   MPI_Comm            comm;
1212259aa2eSHong Zhang   PetscMPIInt         size,rank;
1222259aa2eSHong Zhang   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
1232259aa2eSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1242259aa2eSHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1252259aa2eSHong Zhang 
1262259aa2eSHong Zhang   /* check if matrix local sizes are compatible -- MV! */
1272259aa2eSHong Zhang   if (A->rmap->rstart != P->rmap->rstart || A->rmap->rend != P->rmap->rend) {
1282259aa2eSHong 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);
1292259aa2eSHong Zhang   }
1302259aa2eSHong Zhang   if (A->cmap->rstart != P->rmap->rstart || A->cmap->rend != P->rmap->rend) {
1312259aa2eSHong 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);
1322259aa2eSHong Zhang   }
1332259aa2eSHong Zhang 
13415a3b8e2SHong Zhang   //ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ(A,P,fill,C);CHKERRQ(ierr); //rm later !!!
13515a3b8e2SHong Zhang   //=============================================================================
13615a3b8e2SHong Zhang   Mat                 Cmpi;
13715a3b8e2SHong Zhang   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
138*de817e96SHong Zhang   Mat_SeqAIJ          *p_loc;
139*de817e96SHong Zhang   PetscInt            *pi_loc;
140*de817e96SHong Zhang   PetscInt            *pdti,*pdtj,*poti,*potj,*ptJ,nnz;
14115a3b8e2SHong Zhang   PetscInt            *lnk,*owners_co,*coi,*coj,i,k,pnz,row;
14215a3b8e2SHong Zhang   PetscInt            am=A->rmap->n,pN=P->cmap->N,pm=P->rmap->n,pn=P->cmap->n;
14315a3b8e2SHong Zhang   PetscBT             lnkbt;
14415a3b8e2SHong Zhang   PetscMPIInt         tagi,tagj,*len_si,*len_s,*len_ri,icompleted=0;
14515a3b8e2SHong Zhang   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
14615a3b8e2SHong Zhang   PetscInt            len,proc,*dnz,*onz,*owners;
14715a3b8e2SHong Zhang   PetscInt            nzi,*pti,*ptj;
14815a3b8e2SHong Zhang   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
14915a3b8e2SHong Zhang   MPI_Request         *swaits,*rwaits;
15015a3b8e2SHong Zhang   MPI_Status          *sstatus,rstatus;
15115a3b8e2SHong Zhang   Mat_Merge_SeqsToMPI *merge;
15215a3b8e2SHong Zhang   PetscInt            *api,*apj,*Jptr,apnz,*prmap=p->garray,pon,nspacedouble=0,j,ap_rmax=0;
15315a3b8e2SHong Zhang   PetscReal           afill=1.0,afill_tmp;
15415a3b8e2SHong Zhang   PetscInt            rmax;
15515a3b8e2SHong Zhang 
15615a3b8e2SHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
15715a3b8e2SHong Zhang   ierr        = PetscNew(&ptap);CHKERRQ(ierr);
15815a3b8e2SHong Zhang   ierr        = PetscNew(&merge);CHKERRQ(ierr);
15915a3b8e2SHong Zhang   ptap->merge = merge;
16015a3b8e2SHong Zhang   ptap->reuse = MAT_INITIAL_MATRIX;
16115a3b8e2SHong Zhang 
16215a3b8e2SHong Zhang   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
16315a3b8e2SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
16415a3b8e2SHong Zhang 
16515a3b8e2SHong Zhang   /* get P_loc by taking all local rows of P */
16615a3b8e2SHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
16715a3b8e2SHong Zhang 
168*de817e96SHong Zhang   ptap->reuse = MAT_INITIAL_MATRIX;
16915a3b8e2SHong Zhang 
170*de817e96SHong Zhang   /* (1) compute symbolic AP = A*P, then get AP_loc */
17115a3b8e2SHong Zhang   /*--------------------------------------------------------------------------*/
172*de817e96SHong Zhang   ierr = MatTranspose_SeqAIJ(p->A,MAT_INITIAL_MATRIX,&ptap->Rd);CHKERRQ(ierr);
173*de817e96SHong Zhang   ierr = MatTranspose_SeqAIJ(p->B,MAT_INITIAL_MATRIX,&ptap->Ro);CHKERRQ(ierr);
17415a3b8e2SHong Zhang 
175*de817e96SHong Zhang   ierr = MatMatMult(A,P,MAT_INITIAL_MATRIX,2.0,&AP);CHKERRQ(ierr);
176*de817e96SHong Zhang   ierr = MatMPIAIJGetLocalMat(AP,MAT_INITIAL_MATRIX,&ptap->AP_loc);CHKERRQ(ierr);
177*de817e96SHong Zhang   ierr = MatDestroy(&AP);CHKERRQ(ierr);
17815a3b8e2SHong Zhang 
179*de817e96SHong Zhang   /* (2) compute C_loc=Rd*AP_loc, Co=Ro*AP_loc */
180*de817e96SHong Zhang   ierr = MatMatMult_SeqAIJ_SeqAIJ(ptap->Rd,ptap->AP_loc,MAT_INITIAL_MATRIX,2.0,&ptap->C_loc);CHKERRQ(ierr);
181*de817e96SHong Zhang   ierr = MatMatMult_SeqAIJ_SeqAIJ(ptap->Ro,ptap->AP_loc,MAT_INITIAL_MATRIX,2.0,&ptap->C_oth);CHKERRQ(ierr);
18215a3b8e2SHong Zhang 
18315a3b8e2SHong Zhang 
184*de817e96SHong Zhang   p_loc  = (Mat_SeqAIJ*)(ptap->P_loc)->data;
185*de817e96SHong Zhang   pi_loc = p_loc->i;
18615a3b8e2SHong Zhang 
187*de817e96SHong Zhang   Mat_SeqAIJ *ap = (Mat_SeqAIJ*)ptap->AP_loc->data;
188*de817e96SHong Zhang   api = ap->i;
189*de817e96SHong Zhang   apj = ap->j;
19015a3b8e2SHong Zhang 
19115a3b8e2SHong Zhang   /* (2) determine symbolic Co=(p->B)^T*AP - send to others (coi,coj)*/
19215a3b8e2SHong Zhang   /*-----------------------------------------------------------------*/
19315a3b8e2SHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr);
19415a3b8e2SHong Zhang 
19515a3b8e2SHong Zhang   /* then, compute symbolic Co = (p->B)^T*AP */
19615a3b8e2SHong Zhang   pon    = (p->B)->cmap->n; /* total num of rows to be sent to other processors
19715a3b8e2SHong Zhang                          >= (num of nonzero rows of C_seq) - pn */
19815a3b8e2SHong Zhang   ierr   = PetscMalloc1(pon+1,&coi);CHKERRQ(ierr);
19915a3b8e2SHong Zhang   coi[0] = 0;
20015a3b8e2SHong Zhang 
20115a3b8e2SHong Zhang   /* set initial free space to be fill*(nnz(p->B) + nnz(AP)) */
20215a3b8e2SHong Zhang   nnz           = fill*(poti[pon] + api[am]);
20315a3b8e2SHong Zhang   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
20415a3b8e2SHong Zhang   current_space = free_space;
20515a3b8e2SHong Zhang 
206*de817e96SHong Zhang   /* create and initialize a linked list */
207*de817e96SHong Zhang   ierr = PetscLLCondensedCreate(pN,pN,&lnk,&lnkbt);CHKERRQ(ierr);
208*de817e96SHong Zhang 
20915a3b8e2SHong Zhang   for (i=0; i<pon; i++) {
21015a3b8e2SHong Zhang     pnz = poti[i+1] - poti[i];
21115a3b8e2SHong Zhang     ptJ = potj + poti[i];
21215a3b8e2SHong Zhang     for (j=0; j<pnz; j++) {
21315a3b8e2SHong Zhang       row  = ptJ[j]; /* row of AP == col of Pot */
21415a3b8e2SHong Zhang       apnz = api[row+1] - api[row];
21515a3b8e2SHong Zhang       Jptr = apj + api[row];
21615a3b8e2SHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
21715a3b8e2SHong Zhang       ierr = PetscLLCondensedAddSorted(apnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
21815a3b8e2SHong Zhang     }
21915a3b8e2SHong Zhang     nnz = lnk[0];
22015a3b8e2SHong Zhang 
22115a3b8e2SHong Zhang     /* If free space is not available, double the total space in the list */
22215a3b8e2SHong Zhang     if (current_space->local_remaining<nnz) {
22315a3b8e2SHong Zhang       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
22415a3b8e2SHong Zhang       nspacedouble++;
22515a3b8e2SHong Zhang     }
22615a3b8e2SHong Zhang 
22715a3b8e2SHong Zhang     /* Copy data into free space, and zero out denserows */
22815a3b8e2SHong Zhang     ierr = PetscLLCondensedClean(pN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
22915a3b8e2SHong Zhang 
23015a3b8e2SHong Zhang     current_space->array           += nnz;
23115a3b8e2SHong Zhang     current_space->local_used      += nnz;
23215a3b8e2SHong Zhang     current_space->local_remaining -= nnz;
23315a3b8e2SHong Zhang 
23415a3b8e2SHong Zhang     coi[i+1] = coi[i] + nnz;
23515a3b8e2SHong Zhang   }
23615a3b8e2SHong Zhang 
23715a3b8e2SHong Zhang   ierr      = PetscMalloc1(coi[pon],&coj);CHKERRQ(ierr);
23815a3b8e2SHong Zhang   ierr      = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr);
23915a3b8e2SHong Zhang   afill_tmp = (PetscReal)coi[pon]/(poti[pon] + api[am]+1);
24015a3b8e2SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
24115a3b8e2SHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr);
24215a3b8e2SHong Zhang 
24315a3b8e2SHong Zhang   /* (3) send j-array (coj) of Co to other processors */
24415a3b8e2SHong Zhang   /*--------------------------------------------------*/
24515a3b8e2SHong Zhang   ierr = PetscCalloc1(size,&merge->len_s);CHKERRQ(ierr);
24615a3b8e2SHong Zhang   len_s        = merge->len_s;
24715a3b8e2SHong Zhang   merge->nsend = 0;
24815a3b8e2SHong Zhang 
24915a3b8e2SHong Zhang   /* determine row ownership */
25015a3b8e2SHong Zhang   ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr);
25115a3b8e2SHong Zhang   merge->rowmap->n  = pn;
25215a3b8e2SHong Zhang   merge->rowmap->bs = 1;
25315a3b8e2SHong Zhang 
25415a3b8e2SHong Zhang   ierr   = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr);
25515a3b8e2SHong Zhang   owners = merge->rowmap->range;
25615a3b8e2SHong Zhang 
25715a3b8e2SHong Zhang   /* determine the number of messages to send, their lengths */
25815a3b8e2SHong Zhang   ierr = PetscMalloc2(size,&len_si,size,&sstatus);CHKERRQ(ierr);
25915a3b8e2SHong Zhang   ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
26015a3b8e2SHong Zhang   ierr = PetscMalloc1(size+2,&owners_co);CHKERRQ(ierr);
26115a3b8e2SHong Zhang 
26215a3b8e2SHong Zhang   proc = 0;
26315a3b8e2SHong Zhang   for (i=0; i<pon; i++) {
26415a3b8e2SHong Zhang     while (prmap[i] >= owners[proc+1]) proc++;
26515a3b8e2SHong Zhang     len_si[proc]++;               /* num of rows in Co(=Pt*AP) to be sent to [proc] */
26615a3b8e2SHong Zhang     len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */
26715a3b8e2SHong Zhang   }
26815a3b8e2SHong Zhang 
26915a3b8e2SHong Zhang   len          = 0; /* max length of buf_si[], see (4) */
27015a3b8e2SHong Zhang   owners_co[0] = 0;
27115a3b8e2SHong Zhang   for (proc=0; proc<size; proc++) {
27215a3b8e2SHong Zhang     owners_co[proc+1] = owners_co[proc] + len_si[proc];
27315a3b8e2SHong Zhang     if (len_s[proc]) {
27415a3b8e2SHong Zhang       merge->nsend++;
27515a3b8e2SHong Zhang       len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */
27615a3b8e2SHong Zhang       len         += len_si[proc];
27715a3b8e2SHong Zhang     }
27815a3b8e2SHong Zhang   }
27915a3b8e2SHong Zhang 
28015a3b8e2SHong Zhang   /* determine the number and length of messages to receive for coi and coj  */
28115a3b8e2SHong Zhang   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&merge->nrecv);CHKERRQ(ierr);
28215a3b8e2SHong Zhang   ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr);
28315a3b8e2SHong Zhang 
28415a3b8e2SHong Zhang   /* post the Irecv and Isend of coj */
28515a3b8e2SHong Zhang   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
28615a3b8e2SHong Zhang   ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
28715a3b8e2SHong Zhang   ierr = PetscMalloc1(merge->nsend+1,&swaits);CHKERRQ(ierr);
28815a3b8e2SHong Zhang   for (proc=0, k=0; proc<size; proc++) {
28915a3b8e2SHong Zhang     if (!len_s[proc]) continue;
29015a3b8e2SHong Zhang     i    = owners_co[proc];
29115a3b8e2SHong Zhang     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
29215a3b8e2SHong Zhang     k++;
29315a3b8e2SHong Zhang   }
29415a3b8e2SHong Zhang 
29515a3b8e2SHong Zhang   /* receives and sends of coj are complete */
29615a3b8e2SHong Zhang   for (i=0; i<merge->nrecv; i++) {
29715a3b8e2SHong Zhang     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
29815a3b8e2SHong Zhang   }
29915a3b8e2SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
30015a3b8e2SHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
30115a3b8e2SHong Zhang 
30215a3b8e2SHong Zhang   /* (4) send and recv coi */
30315a3b8e2SHong Zhang   /*-----------------------*/
30415a3b8e2SHong Zhang   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
30515a3b8e2SHong Zhang   ierr   = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
30615a3b8e2SHong Zhang   ierr   = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr);
30715a3b8e2SHong Zhang   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
30815a3b8e2SHong Zhang   for (proc=0,k=0; proc<size; proc++) {
30915a3b8e2SHong Zhang     if (!len_s[proc]) continue;
31015a3b8e2SHong Zhang     /* form outgoing message for i-structure:
31115a3b8e2SHong Zhang          buf_si[0]:                 nrows to be sent
31215a3b8e2SHong Zhang                [1:nrows]:           row index (global)
31315a3b8e2SHong Zhang                [nrows+1:2*nrows+1]: i-structure index
31415a3b8e2SHong Zhang     */
31515a3b8e2SHong Zhang     /*-------------------------------------------*/
31615a3b8e2SHong Zhang     nrows       = len_si[proc]/2 - 1; /* num of rows in Co to be sent to [proc] */
31715a3b8e2SHong Zhang     buf_si_i    = buf_si + nrows+1;
31815a3b8e2SHong Zhang     buf_si[0]   = nrows;
31915a3b8e2SHong Zhang     buf_si_i[0] = 0;
32015a3b8e2SHong Zhang     nrows       = 0;
32115a3b8e2SHong Zhang     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
32215a3b8e2SHong Zhang       nzi = coi[i+1] - coi[i];
32315a3b8e2SHong Zhang       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi;  /* i-structure */
32415a3b8e2SHong Zhang       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
32515a3b8e2SHong Zhang       nrows++;
32615a3b8e2SHong Zhang     }
32715a3b8e2SHong Zhang     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
32815a3b8e2SHong Zhang     k++;
32915a3b8e2SHong Zhang     buf_si += len_si[proc];
33015a3b8e2SHong Zhang   }
33115a3b8e2SHong Zhang   i = merge->nrecv;
33215a3b8e2SHong Zhang   while (i--) {
33315a3b8e2SHong Zhang     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
33415a3b8e2SHong Zhang   }
33515a3b8e2SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
33615a3b8e2SHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
33715a3b8e2SHong Zhang 
33815a3b8e2SHong Zhang   ierr = PetscFree2(len_si,sstatus);CHKERRQ(ierr);
33915a3b8e2SHong Zhang   ierr = PetscFree(len_ri);CHKERRQ(ierr);
34015a3b8e2SHong Zhang   ierr = PetscFree(swaits);CHKERRQ(ierr);
34115a3b8e2SHong Zhang   ierr = PetscFree(buf_s);CHKERRQ(ierr);
34215a3b8e2SHong Zhang 
34315a3b8e2SHong Zhang   /* (5) compute the local portion of C (mpi mat) */
34415a3b8e2SHong Zhang   /*----------------------------------------------*/
34515a3b8e2SHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr);
34615a3b8e2SHong Zhang 
34715a3b8e2SHong Zhang   /* allocate pti array and free space for accumulating nonzero column info */
34815a3b8e2SHong Zhang   ierr   = PetscMalloc1(pn+1,&pti);CHKERRQ(ierr);
34915a3b8e2SHong Zhang   pti[0] = 0;
35015a3b8e2SHong Zhang 
35115a3b8e2SHong Zhang   /* set initial free space to be fill*(nnz(P) + nnz(AP)) */
35215a3b8e2SHong Zhang   nnz           = fill*(pi_loc[pm] + api[am]);
35315a3b8e2SHong Zhang   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
35415a3b8e2SHong Zhang   current_space = free_space;
35515a3b8e2SHong Zhang 
35615a3b8e2SHong Zhang   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
35715a3b8e2SHong Zhang   for (k=0; k<merge->nrecv; k++) {
35815a3b8e2SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
35915a3b8e2SHong Zhang     nrows       = *buf_ri_k[k];
36015a3b8e2SHong Zhang     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
36115a3b8e2SHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
36215a3b8e2SHong Zhang   }
36315a3b8e2SHong Zhang   ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr);
36415a3b8e2SHong Zhang   rmax = 0;
36515a3b8e2SHong Zhang   for (i=0; i<pn; i++) {
36615a3b8e2SHong Zhang     /* add pdt[i,:]*AP into lnk */
36715a3b8e2SHong Zhang     pnz = pdti[i+1] - pdti[i];
36815a3b8e2SHong Zhang     ptJ = pdtj + pdti[i];
36915a3b8e2SHong Zhang     for (j=0; j<pnz; j++) {
37015a3b8e2SHong Zhang       row  = ptJ[j];  /* row of AP == col of Pt */
37115a3b8e2SHong Zhang       apnz = api[row+1] - api[row];
37215a3b8e2SHong Zhang       Jptr = apj + api[row];
37315a3b8e2SHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
37415a3b8e2SHong Zhang       ierr = PetscLLCondensedAddSorted(apnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
37515a3b8e2SHong Zhang     }
37615a3b8e2SHong Zhang 
37715a3b8e2SHong Zhang     /* add received col data into lnk */
37815a3b8e2SHong Zhang     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
37915a3b8e2SHong Zhang       if (i == *nextrow[k]) { /* i-th row */
38015a3b8e2SHong Zhang         nzi  = *(nextci[k]+1) - *nextci[k];
38115a3b8e2SHong Zhang         Jptr = buf_rj[k] + *nextci[k];
38215a3b8e2SHong Zhang         ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr);
38315a3b8e2SHong Zhang         nextrow[k]++; nextci[k]++;
38415a3b8e2SHong Zhang       }
38515a3b8e2SHong Zhang     }
38615a3b8e2SHong Zhang     nnz = lnk[0];
38715a3b8e2SHong Zhang 
38815a3b8e2SHong Zhang     /* if free space is not available, make more free space */
38915a3b8e2SHong Zhang     if (current_space->local_remaining<nnz) {
39015a3b8e2SHong Zhang       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
39115a3b8e2SHong Zhang       nspacedouble++;
39215a3b8e2SHong Zhang     }
39315a3b8e2SHong Zhang     /* copy data into free space, then initialize lnk */
39415a3b8e2SHong Zhang     ierr = PetscLLCondensedClean(pN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
39515a3b8e2SHong Zhang     ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr);
39615a3b8e2SHong Zhang 
39715a3b8e2SHong Zhang     current_space->array           += nnz;
39815a3b8e2SHong Zhang     current_space->local_used      += nnz;
39915a3b8e2SHong Zhang     current_space->local_remaining -= nnz;
40015a3b8e2SHong Zhang 
40115a3b8e2SHong Zhang     pti[i+1] = pti[i] + nnz;
40215a3b8e2SHong Zhang     if (nnz > rmax) rmax = nnz;
40315a3b8e2SHong Zhang   }
40415a3b8e2SHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr);
40515a3b8e2SHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
40615a3b8e2SHong Zhang 
40715a3b8e2SHong Zhang   ierr      = PetscMalloc1(pti[pn]+1,&ptj);CHKERRQ(ierr);
40815a3b8e2SHong Zhang   ierr      = PetscFreeSpaceContiguous(&free_space,ptj);CHKERRQ(ierr);
40915a3b8e2SHong Zhang   afill_tmp = (PetscReal)pti[pn]/(pi_loc[pm] + api[am]+1);
41015a3b8e2SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
41115a3b8e2SHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
41215a3b8e2SHong Zhang 
41315a3b8e2SHong Zhang   /* (6) create symbolic parallel matrix Cmpi */
41415a3b8e2SHong Zhang   /*------------------------------------------*/
41515a3b8e2SHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
41615a3b8e2SHong Zhang   ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
41715a3b8e2SHong Zhang   ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr);
41815a3b8e2SHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
41915a3b8e2SHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
42015a3b8e2SHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
42115a3b8e2SHong Zhang 
42215a3b8e2SHong Zhang   merge->bi        = pti;      /* Cseq->i */
42315a3b8e2SHong Zhang   merge->bj        = ptj;      /* Cseq->j */
42415a3b8e2SHong Zhang   merge->coi       = coi;      /* Co->i   */
42515a3b8e2SHong Zhang   merge->coj       = coj;      /* Co->j   */
42615a3b8e2SHong Zhang   merge->buf_ri    = buf_ri;
42715a3b8e2SHong Zhang   merge->buf_rj    = buf_rj;
42815a3b8e2SHong Zhang   merge->owners_co = owners_co;
42915a3b8e2SHong Zhang   merge->destroy   = Cmpi->ops->destroy;
43015a3b8e2SHong Zhang   merge->duplicate = Cmpi->ops->duplicate;
43115a3b8e2SHong Zhang 
43215a3b8e2SHong Zhang   /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */
43315a3b8e2SHong Zhang   Cmpi->assembled      = PETSC_FALSE;
43415a3b8e2SHong Zhang   Cmpi->ops->destroy   = MatDestroy_MPIAIJ_PtAP;
43515a3b8e2SHong Zhang   Cmpi->ops->duplicate = MatDuplicate_MPIAIJ_MatPtAP;
43615a3b8e2SHong Zhang 
43715a3b8e2SHong Zhang   /* attach the supporting struct to Cmpi for reuse */
43815a3b8e2SHong Zhang   c           = (Mat_MPIAIJ*)Cmpi->data;
43915a3b8e2SHong Zhang   c->ptap     = ptap;
440*de817e96SHong Zhang   ptap->api   = NULL;
441*de817e96SHong Zhang   ptap->apj   = NULL;
44215a3b8e2SHong Zhang   ptap->rmax  = ap_rmax;
44315a3b8e2SHong Zhang   *C          = Cmpi;
44415a3b8e2SHong Zhang 
4452259aa2eSHong Zhang   /* flag 'scalable' determines which implementations to be used:
4462259aa2eSHong Zhang        0: do dense axpy in MatPtAPNumeric() - fast, but requires storage of a nonscalable dense array apa;
4472259aa2eSHong Zhang        1: do sparse axpy in MatPtAPNumeric() - might slow, uses a sparse array apa */
4482259aa2eSHong Zhang   /* set default scalable */
4492259aa2eSHong Zhang   ptap->scalable = PETSC_FALSE; //PETSC_TRUE;
4502259aa2eSHong Zhang 
45115a3b8e2SHong Zhang   ierr = PetscOptionsGetBool(((PetscObject)(*C))->prefix,"-matptap_scalable",&ptap->scalable,NULL);CHKERRQ(ierr);
4522259aa2eSHong Zhang   if (!ptap->scalable) {  /* Do dense axpy */
4532259aa2eSHong Zhang     ierr = PetscCalloc1(P->cmap->N,&ptap->apa);CHKERRQ(ierr);
4542259aa2eSHong Zhang   } else {
4552259aa2eSHong Zhang     //ierr = PetscCalloc1(ap_rmax+1,&ptap->apa);CHKERRQ(ierr);
4562259aa2eSHong Zhang   }
4572259aa2eSHong Zhang 
4580d3441aeSHong Zhang   PetscFunctionReturn(0);
4590d3441aeSHong Zhang }
4600d3441aeSHong Zhang 
4610d3441aeSHong Zhang #undef __FUNCT__
4620d3441aeSHong Zhang #define __FUNCT__ "MatPtAPNumeric_MPIAIJ_MPIAIJ_new"
4630d3441aeSHong Zhang PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ_new(Mat A,Mat P,Mat C)
4640d3441aeSHong Zhang {
4650d3441aeSHong Zhang   PetscErrorCode    ierr;
4660d3441aeSHong Zhang   Mat_MPIAIJ        *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
4670d3441aeSHong Zhang   Mat_SeqAIJ        *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
4680d3441aeSHong Zhang   Mat_PtAPMPI       *ptap = c->ptap;
4699ce11a7cSHong Zhang   Mat               AP_loc,C_loc,C_oth;
4700d3441aeSHong Zhang   PetscInt          i,rstart,rend,cm,ncols,row;
4710d3441aeSHong Zhang   PetscMPIInt       rank;
4720d3441aeSHong Zhang   MPI_Comm          comm;
4730d3441aeSHong Zhang   const PetscInt    *cols;
4740d3441aeSHong Zhang   const PetscScalar *vals;
4750d3441aeSHong Zhang   PetscLogDouble    t0,t1,t2,t3,t4,eR,eAP,eCseq,eCmpi;
4760d3441aeSHong Zhang 
4770d3441aeSHong Zhang   PetscFunctionBegin;
4780d3441aeSHong Zhang   ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr);
4790d3441aeSHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
4800d3441aeSHong Zhang 
4810d3441aeSHong Zhang   ierr = MatZeroEntries(C);CHKERRQ(ierr);
4820d3441aeSHong Zhang 
483e497d3c8SHong Zhang   /* 1) get R = Pd^T,Ro = Po^T */
4840d3441aeSHong Zhang   ierr = PetscTime(&t0);CHKERRQ(ierr);
4850d3441aeSHong Zhang   ierr = MatTranspose_SeqAIJ(p->A,MAT_REUSE_MATRIX,&ptap->Rd);CHKERRQ(ierr);
4860d3441aeSHong Zhang   ierr = MatTranspose_SeqAIJ(p->B,MAT_REUSE_MATRIX,&ptap->Ro);CHKERRQ(ierr);
4870d3441aeSHong Zhang   ierr = PetscTime(&t1);CHKERRQ(ierr);
4880d3441aeSHong Zhang   eR = t1 - t0;
4890d3441aeSHong Zhang 
490e497d3c8SHong Zhang   /* 2) get AP_loc */
4910d3441aeSHong Zhang   AP_loc = ptap->AP_loc;
492e497d3c8SHong Zhang   Mat_SeqAIJ    *ap=(Mat_SeqAIJ*)AP_loc->data;
4930d3441aeSHong Zhang 
494e497d3c8SHong Zhang   /* 2-1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
4950d3441aeSHong Zhang   /*-----------------------------------------------------*/
496e497d3c8SHong Zhang   if (ptap->reuse == MAT_INITIAL_MATRIX) {
4970d3441aeSHong Zhang     /* P_oth and P_loc are obtained in MatPtASymbolic(), skip calling MatGetBrowsOfAoCols() and MatMPIAIJGetLocalMat() */
498e497d3c8SHong Zhang     ptap->reuse = MAT_REUSE_MATRIX;
499e497d3c8SHong Zhang   } else { /* update numerical values of P_oth and P_loc */
5000d3441aeSHong Zhang     ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
5010d3441aeSHong Zhang     ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
502e497d3c8SHong Zhang   }
5030d3441aeSHong Zhang 
504e497d3c8SHong Zhang   /* 2-2) compute numeric C_seq = P_loc^T*A_loc*P - dominating part */
5050d3441aeSHong Zhang   /*--------------------------------------------------------------*/
5060d3441aeSHong Zhang   /* get data from symbolic products */
5070d3441aeSHong Zhang   Mat_SeqAIJ  *p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
5080d3441aeSHong Zhang   Mat_SeqAIJ  *p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
509e497d3c8SHong Zhang   PetscInt    *api,*apj,am = A->rmap->n,j,col,apnz;
510e497d3c8SHong Zhang   PetscScalar *apa = ptap->apa;
5110d3441aeSHong Zhang 
5122259aa2eSHong Zhang   api = ap->i;
5132259aa2eSHong Zhang   apj = ap->j;
514e497d3c8SHong Zhang   for (i=0; i<am; i++) {
515e497d3c8SHong Zhang     /* AP[i,:] = A[i,:]*P = Ad*P_loc Ao*P_oth */
516e497d3c8SHong Zhang     AProw_nonscalable(i,ad,ao,p_loc,p_oth,apa);
517e497d3c8SHong Zhang     apnz = api[i+1] - api[i];
518e497d3c8SHong Zhang     for (j=0; j<apnz; j++) {
519e497d3c8SHong Zhang       col = apj[j+api[i]];
520e497d3c8SHong Zhang       ap->a[j+ap->i[i]] = apa[col];
521e497d3c8SHong Zhang       apa[col] = 0.0;
5220d3441aeSHong Zhang     }
5230d3441aeSHong Zhang   }
5240d3441aeSHong Zhang 
5250d3441aeSHong Zhang   ierr = PetscTime(&t2);CHKERRQ(ierr);
5260d3441aeSHong Zhang   eAP = t2 - t1;
5270d3441aeSHong Zhang 
528e497d3c8SHong Zhang   /* 3) C_loc = R*AP_loc, Co = Ro*AP_loc */
5299ce11a7cSHong Zhang   ierr = MatMatMult_SeqAIJ_SeqAIJ(ptap->Rd,AP_loc,MAT_REUSE_MATRIX,2.0,&ptap->C_loc);CHKERRQ(ierr);
5309ce11a7cSHong Zhang   ierr = MatMatMult_SeqAIJ_SeqAIJ(ptap->Ro,AP_loc,MAT_REUSE_MATRIX,2.0,&ptap->C_oth);CHKERRQ(ierr);
5319ce11a7cSHong Zhang   C_loc = ptap->C_loc;
5329ce11a7cSHong Zhang   C_oth = ptap->C_oth;
5330d3441aeSHong Zhang   //printf("[%d] Co %d, %d\n", rank,Co->rmap->N,Co->cmap->N);
5340d3441aeSHong Zhang   ierr = PetscTime(&t3);CHKERRQ(ierr);
5350d3441aeSHong Zhang   eCseq = t3 - t2;
5360d3441aeSHong Zhang 
5370d3441aeSHong Zhang   /* add C_loc and Co to to C */
5380d3441aeSHong Zhang   ierr = MatGetOwnershipRange(C,&rstart,&rend);CHKERRQ(ierr);
5390d3441aeSHong Zhang 
5400d3441aeSHong Zhang   /* C_loc -> C */
5410d3441aeSHong Zhang   cm = C_loc->rmap->N;
5429ce11a7cSHong Zhang   Mat_SeqAIJ *c_seq;
5439ce11a7cSHong Zhang   c_seq = (Mat_SeqAIJ*)C_loc->data;
5440d3441aeSHong Zhang   for (i=0; i<cm; i++) {
5459ce11a7cSHong Zhang     ncols = c_seq->i[i+1] - c_seq->i[i];
5460d3441aeSHong Zhang     row = rstart + i;
5479ce11a7cSHong Zhang     cols = c_seq->j + c_seq->i[i];
5489ce11a7cSHong Zhang     vals = c_seq->a + c_seq->i[i];
5490d3441aeSHong Zhang     ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr);
5500d3441aeSHong Zhang   }
5510d3441aeSHong Zhang 
5520d3441aeSHong Zhang   /* Co -> C, off-processor part */
5530d3441aeSHong Zhang   //printf("[%d] p->B %d, %d\n",rank,p->B->rmap->N,p->B->cmap->N);
5549ce11a7cSHong Zhang   cm = C_oth->rmap->N;
5559ce11a7cSHong Zhang   c_seq = (Mat_SeqAIJ*)C_oth->data;
5569ce11a7cSHong Zhang   for (i=0; i<cm; i++) {
5579ce11a7cSHong Zhang     ncols = c_seq->i[i+1] - c_seq->i[i];
5580d3441aeSHong Zhang     row = p->garray[i];
5599ce11a7cSHong Zhang     cols = c_seq->j + c_seq->i[i];
5609ce11a7cSHong Zhang     vals = c_seq->a + c_seq->i[i];
5610d3441aeSHong Zhang     //printf("[%d] row[%d] = %d\n",rank,i,row);
5620d3441aeSHong Zhang     ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr);
5630d3441aeSHong Zhang   }
5640d3441aeSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
5650d3441aeSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
5660d3441aeSHong Zhang   ierr = PetscTime(&t4);CHKERRQ(ierr);
5670d3441aeSHong Zhang   eCmpi = t4 - t3;
5680d3441aeSHong Zhang 
5690d3441aeSHong Zhang   if (rank==1) {
5700d3441aeSHong 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);
5710d3441aeSHong Zhang   }
5720d3441aeSHong Zhang   PetscFunctionReturn(0);
5730d3441aeSHong Zhang }
5740d3441aeSHong Zhang 
5750d3441aeSHong Zhang #undef __FUNCT__
57642c57489SHong Zhang #define __FUNCT__ "MatPtAPSymbolic_MPIAIJ_MPIAIJ"
57742c57489SHong Zhang PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
57842c57489SHong Zhang {
57942c57489SHong Zhang   PetscErrorCode      ierr;
58077584fe6SHong Zhang   Mat                 Cmpi;
581f8487c73SHong Zhang   Mat_PtAPMPI         *ptap;
5820298fd71SBarry Smith   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
583f8487c73SHong Zhang   Mat_MPIAIJ          *a        =(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c;
58442c57489SHong Zhang   Mat_SeqAIJ          *ad       =(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
58542c57489SHong Zhang   Mat_SeqAIJ          *p_loc,*p_oth;
586ded4f561SHong Zhang   PetscInt            *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pdti,*pdtj,*poti,*potj,*ptJ;
58777584fe6SHong Zhang   PetscInt            *adi=ad->i,*aj,*aoi=ao->i,nnz;
588a914927fSHong Zhang   PetscInt            *lnk,*owners_co,*coi,*coj,i,k,pnz,row;
589d16ca5f0SHong Zhang   PetscInt            am=A->rmap->n,pN=P->cmap->N,pm=P->rmap->n,pn=P->cmap->n;
59042c57489SHong Zhang   PetscBT             lnkbt;
591ce94432eSBarry Smith   MPI_Comm            comm;
592a914927fSHong Zhang   PetscMPIInt         size,rank,tagi,tagj,*len_si,*len_s,*len_ri,icompleted=0;
59342c57489SHong Zhang   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
59442c57489SHong Zhang   PetscInt            len,proc,*dnz,*onz,*owners;
59524ecddacSHong Zhang   PetscInt            nzi,*pti,*ptj;
59642c57489SHong Zhang   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
597ded4f561SHong Zhang   MPI_Request         *swaits,*rwaits;
598ded4f561SHong Zhang   MPI_Status          *sstatus,rstatus;
59942c57489SHong Zhang   Mat_Merge_SeqsToMPI *merge;
60077584fe6SHong Zhang   PetscInt            *api,*apj,*Jptr,apnz,*prmap=p->garray,pon,nspacedouble=0,j,ap_rmax=0;
601d16ca5f0SHong Zhang   PetscReal           afill=1.0,afill_tmp;
602a914927fSHong Zhang   PetscInt            rmax;
60342c57489SHong Zhang 
60442c57489SHong Zhang   PetscFunctionBegin;
605ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
60624ecddacSHong Zhang 
607c5af039cSHong Zhang   /* check if matrix local sizes are compatible */
608c5af039cSHong Zhang   if (A->rmap->rstart != P->rmap->rstart || A->rmap->rend != P->rmap->rend) {
609c5af039cSHong 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);
610c5af039cSHong Zhang   }
611c5af039cSHong Zhang   if (A->cmap->rstart != P->rmap->rstart || A->cmap->rend != P->rmap->rend) {
612c5af039cSHong 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);
613c5af039cSHong Zhang   }
614c5af039cSHong Zhang 
61583445d95SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
61683445d95SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
61783445d95SHong Zhang 
618f8487c73SHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
619b00a9115SJed Brown   ierr        = PetscNew(&ptap);CHKERRQ(ierr);
620b00a9115SJed Brown   ierr        = PetscNew(&merge);CHKERRQ(ierr);
6219f4155fbSHong Zhang   ptap->merge = merge;
622f8487c73SHong Zhang   ptap->reuse = MAT_INITIAL_MATRIX;
6236c6e00e0SHong Zhang 
6246c6e00e0SHong Zhang   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
625b7f45c76SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
626fe615159SHong Zhang 
6276c6e00e0SHong Zhang   /* get P_loc by taking all local rows of P */
628a1a4e84aSHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
6296c6e00e0SHong Zhang 
630a1a4e84aSHong Zhang   p_loc  = (Mat_SeqAIJ*)(ptap->P_loc)->data;
631a1a4e84aSHong Zhang   p_oth  = (Mat_SeqAIJ*)(ptap->P_oth)->data;
6326c6e00e0SHong Zhang   pi_loc = p_loc->i; pj_loc = p_loc->j;
6336c6e00e0SHong Zhang   pi_oth = p_oth->i; pj_oth = p_oth->j;
63442c57489SHong Zhang 
6352addb229SHong Zhang   /* (1) compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth (api,apj) */
6362addb229SHong Zhang   /*--------------------------------------------------------------------------*/
637854ce69bSBarry Smith   ierr   = PetscMalloc1(am+1,&api);CHKERRQ(ierr);
63882412ba7SHong Zhang   api[0] = 0;
63983445d95SHong Zhang 
64083445d95SHong Zhang   /* create and initialize a linked list */
641a914927fSHong Zhang   ierr = PetscLLCondensedCreate(pN,pN,&lnk,&lnkbt);CHKERRQ(ierr);
64283445d95SHong Zhang 
643a914927fSHong Zhang   /* Initial FreeSpace size is fill*(nnz(A) + nnz(P)) -OOM for ex56, np=8k on Intrepid! */
644d16ca5f0SHong Zhang   ierr = PetscFreeSpaceGet((PetscInt)(fill*(adi[am]+aoi[am]+pi_loc[pm])),&free_space);CHKERRQ(ierr);
6452205254eSKarl Rupp 
646f4a743e1SHong Zhang   current_space = free_space;
647f4a743e1SHong Zhang 
648f4a743e1SHong Zhang   for (i=0; i<am; i++) {
649f4a743e1SHong Zhang     /* diagonal portion of A */
650ded4f561SHong Zhang     nzi = adi[i+1] - adi[i];
65177584fe6SHong Zhang     aj  = ad->j + adi[i];
652ded4f561SHong Zhang     for (j=0; j<nzi; j++) {
65377584fe6SHong Zhang       row  = aj[j];
65482412ba7SHong Zhang       pnz  = pi_loc[row+1] - pi_loc[row];
655ded4f561SHong Zhang       Jptr = pj_loc + pi_loc[row];
65683445d95SHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
657a914927fSHong Zhang       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
658f4a743e1SHong Zhang     }
659f4a743e1SHong Zhang     /* off-diagonal portion of A */
660ded4f561SHong Zhang     nzi = aoi[i+1] - aoi[i];
66177584fe6SHong Zhang     aj  = ao->j + aoi[i];
662ded4f561SHong Zhang     for (j=0; j<nzi; j++) {
66377584fe6SHong Zhang       row  = aj[j];
66482412ba7SHong Zhang       pnz  = pi_oth[row+1] - pi_oth[row];
665ded4f561SHong Zhang       Jptr = pj_oth + pi_oth[row];
666a914927fSHong Zhang       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
667f4a743e1SHong Zhang     }
668a914927fSHong Zhang     apnz     = lnk[0];
66982412ba7SHong Zhang     api[i+1] = api[i] + apnz;
67077584fe6SHong Zhang     if (ap_rmax < apnz) ap_rmax = apnz;
671f4a743e1SHong Zhang 
67283445d95SHong Zhang     /* if free space is not available, double the total space in the list */
67382412ba7SHong Zhang     if (current_space->local_remaining<apnz) {
6742ba1acd5SHong Zhang       ierr = PetscFreeSpaceGet(apnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
675f2b054eeSHong Zhang       nspacedouble++;
676f4a743e1SHong Zhang     }
677f4a743e1SHong Zhang 
678f4a743e1SHong Zhang     /* Copy data into free space, then initialize lnk */
679a914927fSHong Zhang     ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
6802205254eSKarl Rupp 
68182412ba7SHong Zhang     current_space->array           += apnz;
68282412ba7SHong Zhang     current_space->local_used      += apnz;
68382412ba7SHong Zhang     current_space->local_remaining -= apnz;
684f4a743e1SHong Zhang   }
685a914927fSHong Zhang 
68682412ba7SHong Zhang   /* Allocate space for apj, initialize apj, and */
687f4a743e1SHong Zhang   /* destroy list of free space and other temporary array(s) */
688854ce69bSBarry Smith   ierr      = PetscMalloc1(api[am]+1,&apj);CHKERRQ(ierr);
68977584fe6SHong Zhang   ierr      = PetscFreeSpaceContiguous(&free_space,apj);CHKERRQ(ierr);
690118727c9SMark F. Adams   afill_tmp = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1);
691d16ca5f0SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
692f4a743e1SHong Zhang 
6932addb229SHong Zhang   /* (2) determine symbolic Co=(p->B)^T*AP - send to others (coi,coj)*/
6942addb229SHong Zhang   /*-----------------------------------------------------------------*/
695de0260b3SHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr);
69642c57489SHong Zhang 
697ded4f561SHong Zhang   /* then, compute symbolic Co = (p->B)^T*AP */
698d0f46423SBarry Smith   pon    = (p->B)->cmap->n; /* total num of rows to be sent to other processors
69983445d95SHong Zhang                          >= (num of nonzero rows of C_seq) - pn */
700854ce69bSBarry Smith   ierr   = PetscMalloc1(pon+1,&coi);CHKERRQ(ierr);
701de0260b3SHong Zhang   coi[0] = 0;
70242c57489SHong Zhang 
703d16ca5f0SHong Zhang   /* set initial free space to be fill*(nnz(p->B) + nnz(AP)) */
704d16ca5f0SHong Zhang   nnz           = fill*(poti[pon] + api[am]);
70522d28d08SBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
70642c57489SHong Zhang   current_space = free_space;
70742c57489SHong Zhang 
708de0260b3SHong Zhang   for (i=0; i<pon; i++) {
70982412ba7SHong Zhang     pnz = poti[i+1] - poti[i];
71077584fe6SHong Zhang     ptJ = potj + poti[i];
71177584fe6SHong Zhang     for (j=0; j<pnz; j++) {
71277584fe6SHong Zhang       row  = ptJ[j]; /* row of AP == col of Pot */
71382412ba7SHong Zhang       apnz = api[row+1] - api[row];
714ded4f561SHong Zhang       Jptr = apj + api[row];
71583445d95SHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
716a914927fSHong Zhang       ierr = PetscLLCondensedAddSorted(apnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
71742c57489SHong Zhang     }
718a914927fSHong Zhang     nnz = lnk[0];
71942c57489SHong Zhang 
72083445d95SHong Zhang     /* If free space is not available, double the total space in the list */
721ded4f561SHong Zhang     if (current_space->local_remaining<nnz) {
7224238b7adSHong Zhang       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
723d16ca5f0SHong Zhang       nspacedouble++;
72442c57489SHong Zhang     }
72542c57489SHong Zhang 
72642c57489SHong Zhang     /* Copy data into free space, and zero out denserows */
727a914927fSHong Zhang     ierr = PetscLLCondensedClean(pN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
7282205254eSKarl Rupp 
729ded4f561SHong Zhang     current_space->array           += nnz;
730ded4f561SHong Zhang     current_space->local_used      += nnz;
731ded4f561SHong Zhang     current_space->local_remaining -= nnz;
7322205254eSKarl Rupp 
733ded4f561SHong Zhang     coi[i+1] = coi[i] + nnz;
73442c57489SHong Zhang   }
7356b911d16SHong Zhang 
7366b911d16SHong Zhang   ierr      = PetscMalloc1(coi[pon],&coj);CHKERRQ(ierr);
737a1a86e44SBarry Smith   ierr      = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr);
738118727c9SMark F. Adams   afill_tmp = (PetscReal)coi[pon]/(poti[pon] + api[am]+1);
739d16ca5f0SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
740de0260b3SHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr);
74142c57489SHong Zhang 
7426b911d16SHong Zhang   /* (3) send j-array (coj) of Co to other processors */
7436b911d16SHong Zhang   /*--------------------------------------------------*/
7446b911d16SHong Zhang   ierr = PetscCalloc1(size,&merge->len_s);CHKERRQ(ierr);
7456b911d16SHong Zhang   len_s        = merge->len_s;
7466b911d16SHong Zhang   merge->nsend = 0;
7476b911d16SHong Zhang 
7486b911d16SHong Zhang 
74942c57489SHong Zhang   /* determine row ownership */
75026283091SBarry Smith   ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr);
7517a2fc3feSBarry Smith   merge->rowmap->n  = pn;
7527a2fc3feSBarry Smith   merge->rowmap->bs = 1;
7532205254eSKarl Rupp 
75426283091SBarry Smith   ierr   = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr);
7557a2fc3feSBarry Smith   owners = merge->rowmap->range;
75642c57489SHong Zhang 
75742c57489SHong Zhang   /* determine the number of messages to send, their lengths */
758dcca6d9dSJed Brown   ierr = PetscMalloc2(size,&len_si,size,&sstatus);CHKERRQ(ierr);
75983445d95SHong Zhang   ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
760854ce69bSBarry Smith   ierr = PetscMalloc1(size+2,&owners_co);CHKERRQ(ierr);
761de0260b3SHong Zhang 
76283445d95SHong Zhang   proc = 0;
763de0260b3SHong Zhang   for (i=0; i<pon; i++) {
764de0260b3SHong Zhang     while (prmap[i] >= owners[proc+1]) proc++;
7652addb229SHong Zhang     len_si[proc]++;               /* num of rows in Co(=Pt*AP) to be sent to [proc] */
766ee6e4b5bSHong Zhang     len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */
76783445d95SHong Zhang   }
768de0260b3SHong Zhang 
7692addb229SHong Zhang   len          = 0; /* max length of buf_si[], see (4) */
770de0260b3SHong Zhang   owners_co[0] = 0;
77142c57489SHong Zhang   for (proc=0; proc<size; proc++) {
772de0260b3SHong Zhang     owners_co[proc+1] = owners_co[proc] + len_si[proc];
773ee6e4b5bSHong Zhang     if (len_s[proc]) {
77442c57489SHong Zhang       merge->nsend++;
7752addb229SHong Zhang       len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */
77642c57489SHong Zhang       len         += len_si[proc];
77742c57489SHong Zhang     }
77842c57489SHong Zhang   }
77942c57489SHong Zhang 
780ded4f561SHong Zhang   /* determine the number and length of messages to receive for coi and coj  */
7810298fd71SBarry Smith   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&merge->nrecv);CHKERRQ(ierr);
78242c57489SHong Zhang   ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr);
78342c57489SHong Zhang 
784ded4f561SHong Zhang   /* post the Irecv and Isend of coj */
785529bc5dcSHong Zhang   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
786529bc5dcSHong Zhang   ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
787854ce69bSBarry Smith   ierr = PetscMalloc1(merge->nsend+1,&swaits);CHKERRQ(ierr);
78842c57489SHong Zhang   for (proc=0, k=0; proc<size; proc++) {
78942c57489SHong Zhang     if (!len_s[proc]) continue;
790de0260b3SHong Zhang     i    = owners_co[proc];
791529bc5dcSHong Zhang     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
79242c57489SHong Zhang     k++;
79342c57489SHong Zhang   }
79442c57489SHong Zhang 
795ded4f561SHong Zhang   /* receives and sends of coj are complete */
79677584fe6SHong Zhang   for (i=0; i<merge->nrecv; i++) {
797c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
798ded4f561SHong Zhang   }
799ded4f561SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
8000c468ba9SBarry Smith   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
80182412ba7SHong Zhang 
8026b911d16SHong Zhang   /* (4) send and recv coi */
8036b911d16SHong Zhang   /*-----------------------*/
804529bc5dcSHong Zhang   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
805529bc5dcSHong Zhang   ierr   = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
806854ce69bSBarry Smith   ierr   = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr);
80742c57489SHong Zhang   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
80842c57489SHong Zhang   for (proc=0,k=0; proc<size; proc++) {
80942c57489SHong Zhang     if (!len_s[proc]) continue;
81042c57489SHong Zhang     /* form outgoing message for i-structure:
81142c57489SHong Zhang          buf_si[0]:                 nrows to be sent
81242c57489SHong Zhang                [1:nrows]:           row index (global)
81342c57489SHong Zhang                [nrows+1:2*nrows+1]: i-structure index
81442c57489SHong Zhang     */
81542c57489SHong Zhang     /*-------------------------------------------*/
8162addb229SHong Zhang     nrows       = len_si[proc]/2 - 1; /* num of rows in Co to be sent to [proc] */
81742c57489SHong Zhang     buf_si_i    = buf_si + nrows+1;
81842c57489SHong Zhang     buf_si[0]   = nrows;
81942c57489SHong Zhang     buf_si_i[0] = 0;
82042c57489SHong Zhang     nrows       = 0;
821de0260b3SHong Zhang     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
822ded4f561SHong Zhang       nzi = coi[i+1] - coi[i];
823ded4f561SHong Zhang       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi;  /* i-structure */
824de0260b3SHong Zhang       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
82542c57489SHong Zhang       nrows++;
82642c57489SHong Zhang     }
827529bc5dcSHong Zhang     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
82842c57489SHong Zhang     k++;
82942c57489SHong Zhang     buf_si += len_si[proc];
83042c57489SHong Zhang   }
831ded4f561SHong Zhang   i = merge->nrecv;
832ded4f561SHong Zhang   while (i--) {
833c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
834ded4f561SHong Zhang   }
835ded4f561SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
8360c468ba9SBarry Smith   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
837a914927fSHong Zhang 
83824ecddacSHong Zhang   ierr = PetscFree2(len_si,sstatus);CHKERRQ(ierr);
83942c57489SHong Zhang   ierr = PetscFree(len_ri);CHKERRQ(ierr);
840ded4f561SHong Zhang   ierr = PetscFree(swaits);CHKERRQ(ierr);
84142c57489SHong Zhang   ierr = PetscFree(buf_s);CHKERRQ(ierr);
84242c57489SHong Zhang 
8436b911d16SHong Zhang   /* (5) compute the local portion of C (mpi mat) */
8446b911d16SHong Zhang   /*----------------------------------------------*/
845ded4f561SHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr);
846ded4f561SHong Zhang 
84724ecddacSHong Zhang   /* allocate pti array and free space for accumulating nonzero column info */
848854ce69bSBarry Smith   ierr   = PetscMalloc1(pn+1,&pti);CHKERRQ(ierr);
84924ecddacSHong Zhang   pti[0] = 0;
85042c57489SHong Zhang 
851d16ca5f0SHong Zhang   /* set initial free space to be fill*(nnz(P) + nnz(AP)) */
852d16ca5f0SHong Zhang   nnz           = fill*(pi_loc[pm] + api[am]);
85322d28d08SBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
85442c57489SHong Zhang   current_space = free_space;
85542c57489SHong Zhang 
856dcca6d9dSJed Brown   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
85742c57489SHong Zhang   for (k=0; k<merge->nrecv; k++) {
85842c57489SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
85942c57489SHong Zhang     nrows       = *buf_ri_k[k];
86042c57489SHong Zhang     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
86142c57489SHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
86242c57489SHong Zhang   }
86342c57489SHong Zhang   ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr);
86408cb4532SHong Zhang   rmax = 0;
86542c57489SHong Zhang   for (i=0; i<pn; i++) {
866ded4f561SHong Zhang     /* add pdt[i,:]*AP into lnk */
867ded4f561SHong Zhang     pnz = pdti[i+1] - pdti[i];
86877584fe6SHong Zhang     ptJ = pdtj + pdti[i];
86977584fe6SHong Zhang     for (j=0; j<pnz; j++) {
87077584fe6SHong Zhang       row  = ptJ[j];  /* row of AP == col of Pt */
871ded4f561SHong Zhang       apnz = api[row+1] - api[row];
872ded4f561SHong Zhang       Jptr = apj + api[row];
873ded4f561SHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
874a914927fSHong Zhang       ierr = PetscLLCondensedAddSorted(apnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
875ded4f561SHong Zhang     }
876a914927fSHong Zhang 
87742c57489SHong Zhang     /* add received col data into lnk */
87842c57489SHong Zhang     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
87942c57489SHong Zhang       if (i == *nextrow[k]) { /* i-th row */
880ded4f561SHong Zhang         nzi  = *(nextci[k]+1) - *nextci[k];
881ded4f561SHong Zhang         Jptr = buf_rj[k] + *nextci[k];
882a914927fSHong Zhang         ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr);
88342c57489SHong Zhang         nextrow[k]++; nextci[k]++;
88442c57489SHong Zhang       }
88542c57489SHong Zhang     }
886a914927fSHong Zhang     nnz = lnk[0];
88742c57489SHong Zhang 
88842c57489SHong Zhang     /* if free space is not available, make more free space */
889ded4f561SHong Zhang     if (current_space->local_remaining<nnz) {
8904238b7adSHong Zhang       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
891d16ca5f0SHong Zhang       nspacedouble++;
89242c57489SHong Zhang     }
89342c57489SHong Zhang     /* copy data into free space, then initialize lnk */
894a914927fSHong Zhang     ierr = PetscLLCondensedClean(pN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
895ded4f561SHong Zhang     ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr);
8962205254eSKarl Rupp 
897ded4f561SHong Zhang     current_space->array           += nnz;
898ded4f561SHong Zhang     current_space->local_used      += nnz;
899ded4f561SHong Zhang     current_space->local_remaining -= nnz;
9002205254eSKarl Rupp 
90124ecddacSHong Zhang     pti[i+1] = pti[i] + nnz;
90208cb4532SHong Zhang     if (nnz > rmax) rmax = nnz;
90342c57489SHong Zhang   }
904ded4f561SHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr);
9050572522cSBarry Smith   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
90642c57489SHong Zhang 
907854ce69bSBarry Smith   ierr      = PetscMalloc1(pti[pn]+1,&ptj);CHKERRQ(ierr);
90824ecddacSHong Zhang   ierr      = PetscFreeSpaceContiguous(&free_space,ptj);CHKERRQ(ierr);
90924ecddacSHong Zhang   afill_tmp = (PetscReal)pti[pn]/(pi_loc[pm] + api[am]+1);
910d16ca5f0SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
91142c57489SHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
91242c57489SHong Zhang 
9136b911d16SHong Zhang   /* (6) create symbolic parallel matrix Cmpi */
9146b911d16SHong Zhang   /*------------------------------------------*/
91577584fe6SHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
91677584fe6SHong Zhang   ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
91733d57670SJed Brown   ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr);
91877584fe6SHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
91977584fe6SHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
92042c57489SHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
921a2f3521dSMark F. Adams 
922b091e509SHong Zhang   merge->bi        = pti;      /* Cseq->i */
923b091e509SHong Zhang   merge->bj        = ptj;      /* Cseq->j */
924b091e509SHong Zhang   merge->coi       = coi;      /* Co->i   */
925b091e509SHong Zhang   merge->coj       = coj;      /* Co->j   */
92642c57489SHong Zhang   merge->buf_ri    = buf_ri;
92742c57489SHong Zhang   merge->buf_rj    = buf_rj;
928de0260b3SHong Zhang   merge->owners_co = owners_co;
92977584fe6SHong Zhang   merge->destroy   = Cmpi->ops->destroy;
93077584fe6SHong Zhang   merge->duplicate = Cmpi->ops->duplicate;
931cc6cb787SHong Zhang 
93277584fe6SHong Zhang   /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */
933a914927fSHong Zhang   Cmpi->assembled      = PETSC_FALSE;
93477584fe6SHong Zhang   Cmpi->ops->destroy   = MatDestroy_MPIAIJ_PtAP;
93577584fe6SHong Zhang   Cmpi->ops->duplicate = MatDuplicate_MPIAIJ_MatPtAP;
93642c57489SHong Zhang 
93777584fe6SHong Zhang   /* attach the supporting struct to Cmpi for reuse */
93877584fe6SHong Zhang   c           = (Mat_MPIAIJ*)Cmpi->data;
939f8487c73SHong Zhang   c->ptap     = ptap;
94077584fe6SHong Zhang   ptap->api   = api;
94177584fe6SHong Zhang   ptap->apj   = apj;
942d6ab1dc0SHong Zhang   ptap->rmax  = ap_rmax;
94377584fe6SHong Zhang   *C          = Cmpi;
944414894bdSHong Zhang 
945414894bdSHong Zhang   /* flag 'scalable' determines which implementations to be used:
946414894bdSHong Zhang        0: do dense axpy in MatPtAPNumeric() - fast, but requires storage of a nonscalable dense array apa;
947414894bdSHong Zhang        1: do sparse axpy in MatPtAPNumeric() - might slow, uses a sparse array apa */
948414894bdSHong Zhang   /* set default scalable */
9499516a85cSHong Zhang   ptap->scalable = PETSC_FALSE; //PETSC_TRUE;
9502205254eSKarl Rupp 
9510298fd71SBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)Cmpi)->prefix,"-matptap_scalable",&ptap->scalable,NULL);CHKERRQ(ierr);
952414894bdSHong Zhang   if (!ptap->scalable) {  /* Do dense axpy */
9531795a4d1SJed Brown     ierr = PetscCalloc1(pN,&ptap->apa);CHKERRQ(ierr);
954414894bdSHong Zhang   } else {
9551795a4d1SJed Brown     ierr = PetscCalloc1(ap_rmax+1,&ptap->apa);CHKERRQ(ierr);
956414894bdSHong Zhang   }
957414894bdSHong Zhang 
958f2b054eeSHong Zhang #if defined(PETSC_USE_INFO)
95924ecddacSHong Zhang   if (pti[pn] != 0) {
96057622a8eSBarry Smith     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
96157622a8eSBarry Smith     ierr = PetscInfo1(Cmpi,"Use MatPtAP(A,P,MatReuse,%g,&C) for best performance.\n",(double)afill);CHKERRQ(ierr);
962f2b054eeSHong Zhang   } else {
96377584fe6SHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
964f2b054eeSHong Zhang   }
965f2b054eeSHong Zhang #endif
96642c57489SHong Zhang   PetscFunctionReturn(0);
96742c57489SHong Zhang }
96842c57489SHong Zhang 
96942c57489SHong Zhang #undef __FUNCT__
97042c57489SHong Zhang #define __FUNCT__ "MatPtAPNumeric_MPIAIJ_MPIAIJ"
97142c57489SHong Zhang PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C)
97242c57489SHong Zhang {
97342c57489SHong Zhang   PetscErrorCode      ierr;
974f8487c73SHong Zhang   Mat_MPIAIJ          *a =(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
97542c57489SHong Zhang   Mat_SeqAIJ          *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
976de0260b3SHong Zhang   Mat_SeqAIJ          *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data;
97742c57489SHong Zhang   Mat_SeqAIJ          *p_loc,*p_oth;
978f8487c73SHong Zhang   Mat_PtAPMPI         *ptap;
9799f4155fbSHong Zhang   Mat_Merge_SeqsToMPI *merge;
9801d633602SHong Zhang   PetscInt            *adi=ad->i,*aoi=ao->i,*adj,*aoj,*apJ,nextp;
98182412ba7SHong Zhang   PetscInt            *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pJ,*pj;
98282412ba7SHong Zhang   PetscInt            i,j,k,anz,pnz,apnz,nextap,row,*cj;
983e900a757SHong Zhang   MatScalar           *ada,*aoa,*apa,*pa,*ca,*pa_loc,*pa_oth,valtmp;
984d0f46423SBarry Smith   PetscInt            am  =A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n;
985ce94432eSBarry Smith   MPI_Comm            comm;
98642c57489SHong Zhang   PetscMPIInt         size,rank,taga,*len_s;
98782412ba7SHong Zhang   PetscInt            *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci;
98842c57489SHong Zhang   PetscInt            **buf_ri,**buf_rj;
98950f5bf86SHong Zhang   PetscInt            cnz=0,*bj_i,*bi,*bj,bnz,nextcj;  /* bi,bj,ba: local array of C(mpi mat) */
99042c57489SHong Zhang   MPI_Request         *s_waits,*r_waits;
99142c57489SHong Zhang   MPI_Status          *status;
99282412ba7SHong Zhang   MatScalar           **abuf_r,*ba_i,*pA,*coa,*ba;
99382412ba7SHong Zhang   PetscInt            *api,*apj,*coi,*coj;
994d0f46423SBarry Smith   PetscInt            *poJ=po->j,*pdJ=pd->j,pcstart=P->cmap->rstart,pcend=P->cmap->rend;
9956ce94e8fSHong Zhang   PetscBool           scalable;
99638571addSHong Zhang #if defined(PTAP_PROFILE)
997e497d3c8SHong Zhang   PetscLogDouble t0,t1,t2,eP,t3,t4,et2_AP=0.0,ePtAP=0.0,t2_0,t2_1,t2_2;
99838571addSHong Zhang #endif
99942c57489SHong Zhang 
100042c57489SHong Zhang   PetscFunctionBegin;
1001ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr);
100238571addSHong Zhang #if defined(PTAP_PROFILE)
10038563dfccSBarry Smith   ierr = PetscTime(&t0);CHKERRQ(ierr);
100438571addSHong Zhang #endif
100542c57489SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
100642c57489SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
100742c57489SHong Zhang 
1008f8487c73SHong Zhang   ptap = c->ptap;
1009ce94432eSBarry 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");
1010f8487c73SHong Zhang   merge    = ptap->merge;
1011414894bdSHong Zhang   apa      = ptap->apa;
10126ce94e8fSHong Zhang   scalable = ptap->scalable;
10136c6e00e0SHong Zhang 
1014a1a4e84aSHong Zhang   /* 1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
10156b911d16SHong Zhang   /*-----------------------------------------------------*/
1016f8487c73SHong Zhang   if (ptap->reuse == MAT_INITIAL_MATRIX) {
10179f4155fbSHong Zhang     /* P_oth and P_loc are obtained in MatPtASymbolic(), skip calling MatGetBrowsOfAoCols() and MatMPIAIJGetLocalMat() */
1018f8487c73SHong Zhang     ptap->reuse = MAT_REUSE_MATRIX;
10196c6e00e0SHong Zhang   } else { /* update numerical values of P_oth and P_loc */
1020b7f45c76SHong Zhang     ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
1021a1a4e84aSHong Zhang     ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
10226c6e00e0SHong Zhang   }
102338571addSHong Zhang #if defined(PTAP_PROFILE)
10248563dfccSBarry Smith   ierr = PetscTime(&t1);CHKERRQ(ierr);
102505d62848SHong Zhang   eP = t1-t0;
102638571addSHong Zhang #endif
102705d62848SHong Zhang   /*
102805d62848SHong Zhang   printf("[%d] Ad: %d, %d; Ao: %d, %d; P_loc: %d, %d; P_oth %d, %d;\n",rank,
102905d62848SHong Zhang          a->A->rmap->N,a->A->cmap->N,a->B->rmap->N,a->B->cmap->N,
103005d62848SHong Zhang          ptap->P_loc->rmap->N,ptap->P_loc->cmap->N,
103105d62848SHong Zhang          ptap->P_oth->rmap->N,ptap->P_oth->cmap->N);
103205d62848SHong Zhang    */
1033f8487c73SHong Zhang 
1034f9473cf7SHong Zhang   /* 2) compute numeric C_seq = P_loc^T*A_loc*P - dominating part */
1035f9473cf7SHong Zhang   /*--------------------------------------------------------------*/
103642c57489SHong Zhang   /* get data from symbolic products */
1037a1a4e84aSHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
1038a1a4e84aSHong Zhang   p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
1039a9d06231SHong Zhang   pi_loc=p_loc->i; pj_loc=p_loc->j; pJ=pj_loc; pa_loc=p_loc->a;
104042c57489SHong Zhang   pi_oth=p_oth->i; pj_oth=p_oth->j; pa_oth=p_oth->a;
104142c57489SHong Zhang 
1042de0260b3SHong Zhang   coi  = merge->coi; coj = merge->coj;
10431795a4d1SJed Brown   ierr = PetscCalloc1(coi[pon]+1,&coa);CHKERRQ(ierr);
1044de0260b3SHong Zhang 
1045de0260b3SHong Zhang   bi     = merge->bi; bj = merge->bj;
10467a2fc3feSBarry Smith   owners = merge->rowmap->range;
10471795a4d1SJed Brown   ierr   = PetscCalloc1(bi[cm]+1,&ba);CHKERRQ(ierr);  /* ba: Cseq->a */
104842c57489SHong Zhang 
1049a1a4e84aSHong Zhang   api = ptap->api; apj = ptap->apj;
1050d9824c15SHong Zhang 
10519516a85cSHong Zhang   if (!scalable) { /* Do dense axpy on apa (length of pN, stores A[i,:]*P) - nonscalable, but faster (could take 1/3 scalable time) */
1052b38d1a2dSJed Brown     ierr = PetscInfo(C,"Using non-scalable dense axpy\n");CHKERRQ(ierr);
105305d62848SHong Zhang #if 0
10540d3441aeSHong Zhang     /* ------ 10x slower -------------- */
10550d3441aeSHong Zhang     /*==================================*/
10569516a85cSHong Zhang     Mat         R = ptap->R;
10579516a85cSHong Zhang     Mat_SeqAIJ  *r = (Mat_SeqAIJ*)R->data;
10589516a85cSHong Zhang     PetscInt    *ri=r->i,*rj=r->j,rnz,arow,l,prow,pcol,pN=P->cmap->N;
10599516a85cSHong Zhang     PetscScalar *ra=r->a,tmp,cdense[pN];
10609516a85cSHong Zhang 
10619516a85cSHong Zhang     ierr = PetscMemzero(cdense,pN*sizeof(PetscScalar));CHKERRQ(ierr);
10629516a85cSHong Zhang     for (i=0; i<cm; i++) { /* each row of C or R */
10639516a85cSHong Zhang       rnz = ri[i+1] - ri[i];
10649516a85cSHong Zhang 
10659516a85cSHong Zhang       for (j=0; j<rnz; j++) { /* each nz of R */
10669516a85cSHong Zhang         arow = rj[ri[i] + j];
10679516a85cSHong Zhang 
10689516a85cSHong Zhang         /* diagonal portion of A */
10699516a85cSHong Zhang         anz  = ad->i[arow+1] - ad->i[arow];
10709516a85cSHong Zhang         for (k=0; k<anz; k++) { /* each nz of Ad */
10719516a85cSHong Zhang           tmp  = ra[ri[i] + j]*ad->a[ad->i[arow] + k];
10729516a85cSHong Zhang           prow = ad->j[ad->i[arow] + k];
10739516a85cSHong Zhang           pnz  = pi_loc[prow+1] - pi_loc[prow];
10749516a85cSHong Zhang 
10759516a85cSHong Zhang           for (l=0; l<pnz; l++) { /* each nz of P_loc */
10769516a85cSHong Zhang             pcol = pj_loc[pi_loc[prow] + l];
10779516a85cSHong Zhang             cdense[pcol] += tmp*pa_loc[pi_loc[prow] + l];
10789516a85cSHong Zhang           }
10799516a85cSHong Zhang         }
10809516a85cSHong Zhang 
10819516a85cSHong Zhang         /* off-diagonal portion of A */
10829516a85cSHong Zhang         anz  = ao->i[arow+1] - ao->i[arow];
10839516a85cSHong Zhang         for (k=0; k<anz; k++) { /* each nz of Ao */
10849516a85cSHong Zhang           tmp  = ra[ri[i] + j]*ao->a[ao->i[arow] + k];
10859516a85cSHong Zhang           prow = ao->j[ao->i[arow] + k];
10869516a85cSHong Zhang           pnz  = pi_oth[prow+1] - pi_oth[prow];
10879516a85cSHong Zhang 
10889516a85cSHong Zhang           for (l=0; l<pnz; l++) { /* each nz of P_oth */
10899516a85cSHong Zhang             pcol = pj_oth[pi_oth[prow] + l];
10909516a85cSHong Zhang             cdense[pcol] += tmp*pa_oth[pi_oth[prow] + l];
10919516a85cSHong Zhang           }
10929516a85cSHong Zhang         }
10939516a85cSHong Zhang 
10949516a85cSHong Zhang       } //for (j=0; j<rnz; j++)
10959516a85cSHong Zhang 
10969516a85cSHong Zhang       /* copy cdense[] into ca; zero cdense[] */
10979516a85cSHong Zhang       cnz = bi[i+1] - bi[i];
10989516a85cSHong Zhang       cj  = bj + bi[i];
10999516a85cSHong Zhang       ca  = ba + bi[i];
11009516a85cSHong Zhang       for (j=0; j<cnz; j++) {
11019516a85cSHong Zhang         ca[j] += cdense[cj[j]];
11029516a85cSHong Zhang         cdense[cj[j]] = 0.0;
11039516a85cSHong Zhang       }
11049516a85cSHong Zhang #if 0
11059516a85cSHong Zhang       if (rank == 0) {
11069516a85cSHong Zhang         printf("[%d] row %d: ",rank,i);
11079516a85cSHong Zhang         for (j=0; j<pN; j++) printf(" %g,",cdense[j]);
11089516a85cSHong Zhang         printf("\n");
11099516a85cSHong Zhang       }
11109516a85cSHong Zhang       for (j=0; j<pN; j++) cdense[j]=0.0; // zero cdnese[]
11119516a85cSHong Zhang #endif
11129516a85cSHong Zhang     } //for (i=0; i<cm; i++) {
111305d62848SHong Zhang #endif
11149516a85cSHong Zhang 
11159516a85cSHong Zhang     //==========================================
11169516a85cSHong Zhang 
111705d62848SHong Zhang     ierr = PetscTime(&t1);CHKERRQ(ierr);
1118a9d06231SHong Zhang     for (i=0; i<am; i++) {
1119b091e509SHong Zhang #if defined(PTAP_PROFILE)
11208563dfccSBarry Smith       ierr = PetscTime(&t2_0);CHKERRQ(ierr);
1121b091e509SHong Zhang #endif
1122a9d06231SHong Zhang       /* 2-a) form i-th sparse row of A_loc*P = Ad*P_loc + Ao*P_oth */
1123a9d06231SHong Zhang       /*------------------------------------------------------------*/
1124a9d06231SHong Zhang       apJ = apj + api[i];
1125a9d06231SHong Zhang 
1126a9d06231SHong Zhang       /* diagonal portion of A */
1127a9d06231SHong Zhang       anz = adi[i+1] - adi[i];
1128a9d06231SHong Zhang       adj = ad->j + adi[i];
1129a9d06231SHong Zhang       ada = ad->a + adi[i];
1130a9d06231SHong Zhang       for (j=0; j<anz; j++) {
1131a9d06231SHong Zhang         row = adj[j];
1132a9d06231SHong Zhang         pnz = pi_loc[row+1] - pi_loc[row];
1133a9d06231SHong Zhang         pj  = pj_loc + pi_loc[row];
1134a9d06231SHong Zhang         pa  = pa_loc + pi_loc[row];
1135a9d06231SHong Zhang 
1136a9d06231SHong Zhang         /* perform dense axpy */
1137e900a757SHong Zhang         valtmp = ada[j];
1138a9d06231SHong Zhang         for (k=0; k<pnz; k++) {
1139e900a757SHong Zhang           apa[pj[k]] += valtmp*pa[k];
1140a9d06231SHong Zhang         }
1141a9d06231SHong Zhang         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
1142a9d06231SHong Zhang       }
1143a9d06231SHong Zhang 
1144a9d06231SHong Zhang       /* off-diagonal portion of A */
1145a9d06231SHong Zhang       anz = aoi[i+1] - aoi[i];
1146a9d06231SHong Zhang       aoj = ao->j + aoi[i];
1147a9d06231SHong Zhang       aoa = ao->a + aoi[i];
1148a9d06231SHong Zhang       for (j=0; j<anz; j++) {
1149a9d06231SHong Zhang         row = aoj[j];
1150a9d06231SHong Zhang         pnz = pi_oth[row+1] - pi_oth[row];
1151a9d06231SHong Zhang         pj  = pj_oth + pi_oth[row];
1152a9d06231SHong Zhang         pa  = pa_oth + pi_oth[row];
1153a9d06231SHong Zhang 
1154a9d06231SHong Zhang         /* perform dense axpy */
1155e900a757SHong Zhang         valtmp = aoa[j];
1156a9d06231SHong Zhang         for (k=0; k<pnz; k++) {
1157e900a757SHong Zhang           apa[pj[k]] += valtmp*pa[k];
1158a9d06231SHong Zhang         }
1159a9d06231SHong Zhang         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
1160a9d06231SHong Zhang       }
1161b091e509SHong Zhang #if defined(PTAP_PROFILE)
11628563dfccSBarry Smith       ierr    = PetscTime(&t2_1);CHKERRQ(ierr);
1163b091e509SHong Zhang       et2_AP += t2_1 - t2_0;
1164b091e509SHong Zhang #endif
1165a9d06231SHong Zhang 
1166a9d06231SHong Zhang       /* 2-b) Compute Cseq = P_loc[i,:]^T*AP[i,:] using outer product */
1167a9d06231SHong Zhang       /*--------------------------------------------------------------*/
1168a9d06231SHong Zhang       apnz = api[i+1] - api[i];
1169a9d06231SHong Zhang       /* put the value into Co=(p->B)^T*AP (off-diagonal part, send to others) */
1170a9d06231SHong Zhang       pnz = po->i[i+1] - po->i[i];
1171e900a757SHong Zhang       poJ = po->j + po->i[i];
1172a9d06231SHong Zhang       pA  = po->a + po->i[i];
1173a9d06231SHong Zhang       for (j=0; j<pnz; j++) {
1174e900a757SHong Zhang         row = poJ[j];
1175e900a757SHong Zhang         cnz = coi[row+1] - coi[row];
1176e900a757SHong Zhang         cj  = coj + coi[row];
1177e900a757SHong Zhang         ca  = coa + coi[row];
1178a9d06231SHong Zhang         /* perform dense axpy */
1179e900a757SHong Zhang         valtmp = pA[j];
1180a9d06231SHong Zhang         for (k=0; k<cnz; k++) {
1181e900a757SHong Zhang           ca[k] += valtmp*apa[cj[k]];
1182a9d06231SHong Zhang         }
1183a9d06231SHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1184a9d06231SHong Zhang       }
118505d62848SHong Zhang #if 1
1186a9d06231SHong Zhang       /* put the value into Cd (diagonal part) */
1187a9d06231SHong Zhang       pnz = pd->i[i+1] - pd->i[i];
1188e900a757SHong Zhang       pdJ = pd->j + pd->i[i];
1189a9d06231SHong Zhang       pA  = pd->a + pd->i[i];
1190a9d06231SHong Zhang       for (j=0; j<pnz; j++) {
1191e900a757SHong Zhang         row = pdJ[j];
1192e900a757SHong Zhang         cnz = bi[row+1] - bi[row];
1193e900a757SHong Zhang         cj  = bj + bi[row];
1194e900a757SHong Zhang         ca  = ba + bi[row];
1195a9d06231SHong Zhang         /* perform dense axpy */
1196e900a757SHong Zhang         valtmp = pA[j];
1197a9d06231SHong Zhang         for (k=0; k<cnz; k++) {
1198e900a757SHong Zhang           ca[k] += valtmp*apa[cj[k]];
1199a9d06231SHong Zhang         }
1200a9d06231SHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1201a9d06231SHong Zhang       }
12029516a85cSHong Zhang #endif
1203a9d06231SHong Zhang       /* zero the current row of A*P */
1204a9d06231SHong Zhang       for (k=0; k<apnz; k++) apa[apJ[k]] = 0.0;
1205b091e509SHong Zhang #if defined(PTAP_PROFILE)
12068563dfccSBarry Smith       ierr      = PetscTime(&t2_2);CHKERRQ(ierr);
120705d62848SHong Zhang       ePtAP += t2_2 - t2_1;
1208b091e509SHong Zhang #endif
1209a9d06231SHong Zhang     }
12109516a85cSHong Zhang 
12119516a85cSHong Zhang     if (rank == 100) {
12129516a85cSHong Zhang     for (row=0; row<cm; row++) {
12139516a85cSHong Zhang       printf("[%d] row %d: ",rank,row);
12149516a85cSHong Zhang       cnz = bi[row+1] - bi[row];
12159516a85cSHong Zhang       for (j=0; j<cnz; j++) printf(" %g,",ba[bi[row]+j]);
12169516a85cSHong Zhang       printf("\n");
12179516a85cSHong Zhang     }
12189516a85cSHong Zhang     }
12199516a85cSHong Zhang 
122038571addSHong Zhang   } else { /* Do sparse axpy on apa (length of ap_rmax, stores A[i,:]*P) - scalable, but slower */
1221b38d1a2dSJed Brown     ierr = PetscInfo(C,"Using scalable sparse axpy\n");CHKERRQ(ierr);
122238571addSHong Zhang     /*-----------------------------------------------------------------------------------------*/
1223a9d06231SHong Zhang     pA=pa_loc;
122442c57489SHong Zhang     for (i=0; i<am; i++) {
1225b091e509SHong Zhang #if defined(PTAP_PROFILE)
12268563dfccSBarry Smith       ierr = PetscTime(&t2_0);CHKERRQ(ierr);
1227b091e509SHong Zhang #endif
1228f4a743e1SHong Zhang       /* form i-th sparse row of A*P */
122982412ba7SHong Zhang       apnz = api[i+1] - api[i];
123082412ba7SHong Zhang       apJ  = apj + api[i];
123142c57489SHong Zhang       /* diagonal portion of A */
123282412ba7SHong Zhang       anz = adi[i+1] - adi[i];
12331d633602SHong Zhang       adj = ad->j + adi[i];
12341d633602SHong Zhang       ada = ad->a + adi[i];
123582412ba7SHong Zhang       for (j=0; j<anz; j++) {
12361d633602SHong Zhang         row    = adj[j];
123782412ba7SHong Zhang         pnz    = pi_loc[row+1] - pi_loc[row];
123882412ba7SHong Zhang         pj     = pj_loc + pi_loc[row];
123982412ba7SHong Zhang         pa     = pa_loc + pi_loc[row];
1240e900a757SHong Zhang         valtmp = ada[j];
1241f4a743e1SHong Zhang         nextp  = 0;
124282412ba7SHong Zhang         for (k=0; nextp<pnz; k++) {
124382412ba7SHong Zhang           if (apJ[k] == pj[nextp]) { /* col of AP == col of P */
1244e900a757SHong Zhang             apa[k] += valtmp*pa[nextp++];
124542c57489SHong Zhang           }
124642c57489SHong Zhang         }
1247dc0b31edSSatish Balay         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
124842c57489SHong Zhang       }
124942c57489SHong Zhang       /* off-diagonal portion of A */
125082412ba7SHong Zhang       anz = aoi[i+1] - aoi[i];
12511d633602SHong Zhang       aoj = ao->j + aoi[i];
12521d633602SHong Zhang       aoa = ao->a + aoi[i];
125382412ba7SHong Zhang       for (j=0; j<anz; j++) {
12541d633602SHong Zhang         row    = aoj[j];
125582412ba7SHong Zhang         pnz    = pi_oth[row+1] - pi_oth[row];
125682412ba7SHong Zhang         pj     = pj_oth + pi_oth[row];
125782412ba7SHong Zhang         pa     = pa_oth + pi_oth[row];
1258e900a757SHong Zhang         valtmp = aoa[j];
1259f4a743e1SHong Zhang         nextp  = 0;
126082412ba7SHong Zhang         for (k=0; nextp<pnz; k++) {
126182412ba7SHong Zhang           if (apJ[k] == pj[nextp]) { /* col of AP == col of P */
1262e900a757SHong Zhang             apa[k] += valtmp*pa[nextp++];
126342c57489SHong Zhang           }
126442c57489SHong Zhang         }
1265dc0b31edSSatish Balay         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
126642c57489SHong Zhang       }
1267b091e509SHong Zhang #if defined(PTAP_PROFILE)
12688563dfccSBarry Smith       ierr    = PetscTime(&t2_1);CHKERRQ(ierr);
1269b091e509SHong Zhang       et2_AP += t2_1 - t2_0;
1270b091e509SHong Zhang #endif
127142c57489SHong Zhang 
1272a9d06231SHong Zhang       /* 2-b) Compute Cseq = P_loc[i,:]^T*AP[i,:] using outer product */
1273a9d06231SHong Zhang       /*--------------------------------------------------------------*/
127482412ba7SHong Zhang       pnz = pi_loc[i+1] - pi_loc[i];
1275f9473cf7SHong Zhang       pJ  = pj_loc + pi_loc[i];
127682412ba7SHong Zhang       for (j=0; j<pnz; j++) {
127742c57489SHong Zhang         nextap = 0;
1278f9473cf7SHong Zhang         row    = pJ[j]; /* global index */
127982412ba7SHong Zhang         if (row < pcstart || row >=pcend) { /* put the value into Co */
1280e900a757SHong Zhang           row = *poJ;
1281e900a757SHong Zhang           cj  = coj + coi[row];
1282e900a757SHong Zhang           ca  = coa + coi[row]; poJ++;
128382412ba7SHong Zhang         } else {                            /* put the value into Cd */
1284e900a757SHong Zhang           row = *pdJ;
1285e900a757SHong Zhang           cj  = bj + bi[row];
1286e900a757SHong Zhang           ca  = ba + bi[row]; pdJ++;
128742c57489SHong Zhang         }
1288e900a757SHong Zhang         valtmp = pA[j];
128982412ba7SHong Zhang         for (k=0; nextap<apnz; k++) {
1290e900a757SHong Zhang           if (cj[k]==apJ[nextap]) ca[k] += valtmp*apa[nextap++];
129142c57489SHong Zhang         }
1292dc0b31edSSatish Balay         ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr);
1293de0260b3SHong Zhang       }
12940496f32cSHong Zhang       pA += pnz;
129542c57489SHong Zhang       /* zero the current row info for A*P */
129682412ba7SHong Zhang       ierr = PetscMemzero(apa,apnz*sizeof(MatScalar));CHKERRQ(ierr);
1297b091e509SHong Zhang #if defined(PTAP_PROFILE)
12988563dfccSBarry Smith       ierr      = PetscTime(&t2_2);CHKERRQ(ierr);
129905d62848SHong Zhang       ePtAP += t2_2 - t2_1;
1300b091e509SHong Zhang #endif
130142c57489SHong Zhang     }
1302d9824c15SHong Zhang   }
130338571addSHong Zhang #if defined(PTAP_PROFILE)
13048563dfccSBarry Smith   ierr = PetscTime(&t2);CHKERRQ(ierr);
130538571addSHong Zhang #endif
130642c57489SHong Zhang 
1307a9d06231SHong Zhang   /* 3) send and recv matrix values coa */
1308a9d06231SHong Zhang   /*------------------------------------*/
130942c57489SHong Zhang   buf_ri = merge->buf_ri;
131042c57489SHong Zhang   buf_rj = merge->buf_rj;
131142c57489SHong Zhang   len_s  = merge->len_s;
1312fc42d0c8SSatish Balay   ierr   = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr);
131342c57489SHong Zhang   ierr   = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr);
131442c57489SHong Zhang 
1315dcca6d9dSJed Brown   ierr = PetscMalloc2(merge->nsend+1,&s_waits,size,&status);CHKERRQ(ierr);
131642c57489SHong Zhang   for (proc=0,k=0; proc<size; proc++) {
131742c57489SHong Zhang     if (!len_s[proc]) continue;
1318de0260b3SHong Zhang     i    = merge->owners_co[proc];
1319de0260b3SHong Zhang     ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr);
132042c57489SHong Zhang     k++;
132142c57489SHong Zhang   }
13220c468ba9SBarry Smith   if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);}
13230c468ba9SBarry Smith   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);}
132442c57489SHong Zhang 
1325a9d06231SHong Zhang   ierr = PetscFree2(s_waits,status);CHKERRQ(ierr);
132642c57489SHong Zhang   ierr = PetscFree(r_waits);CHKERRQ(ierr);
132782412ba7SHong Zhang   ierr = PetscFree(coa);CHKERRQ(ierr);
132838571addSHong Zhang #if defined(PTAP_PROFILE)
13298563dfccSBarry Smith   ierr = PetscTime(&t3);CHKERRQ(ierr);
133038571addSHong Zhang #endif
133142c57489SHong Zhang 
1332a9d06231SHong Zhang   /* 4) insert local Cseq and received values into Cmpi */
1333a9d06231SHong Zhang   /*------------------------------------------------------*/
1334dcca6d9dSJed Brown   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
133542c57489SHong Zhang   for (k=0; k<merge->nrecv; k++) {
133642c57489SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
133742c57489SHong Zhang     nrows       = *(buf_ri_k[k]);
133842c57489SHong Zhang     nextrow[k]  = buf_ri_k[k]+1;  /* next row number of k-th recved i-structure */
133982412ba7SHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
134042c57489SHong Zhang   }
134142c57489SHong Zhang 
1342de0260b3SHong Zhang   for (i=0; i<cm; i++) {
134382412ba7SHong Zhang     row  = owners[rank] + i; /* global row index of C_seq */
134442c57489SHong Zhang     bj_i = bj + bi[i];  /* col indices of the i-th row of C */
1345de0260b3SHong Zhang     ba_i = ba + bi[i];
134682412ba7SHong Zhang     bnz  = bi[i+1] - bi[i];
134742c57489SHong Zhang     /* add received vals into ba */
134842c57489SHong Zhang     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
134942c57489SHong Zhang       /* i-th row */
135042c57489SHong Zhang       if (i == *nextrow[k]) {
135182412ba7SHong Zhang         cnz    = *(nextci[k]+1) - *nextci[k];
135282412ba7SHong Zhang         cj     = buf_rj[k] + *(nextci[k]);
135382412ba7SHong Zhang         ca     = abuf_r[k] + *(nextci[k]);
135482412ba7SHong Zhang         nextcj = 0;
135582412ba7SHong Zhang         for (j=0; nextcj<cnz; j++) {
135682412ba7SHong Zhang           if (bj_i[j] == cj[nextcj]) { /* bcol == ccol */
135782412ba7SHong Zhang             ba_i[j] += ca[nextcj++];
135842c57489SHong Zhang           }
135942c57489SHong Zhang         }
136082412ba7SHong Zhang         nextrow[k]++; nextci[k]++;
1361c5af039cSHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
136242c57489SHong Zhang       }
136342c57489SHong Zhang     }
136482412ba7SHong Zhang     ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr);
136542c57489SHong Zhang   }
136642c57489SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
136742c57489SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
136842c57489SHong Zhang 
1369de0260b3SHong Zhang   ierr = PetscFree(ba);CHKERRQ(ierr);
1370c05d87d6SBarry Smith   ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr);
137142c57489SHong Zhang   ierr = PetscFree(abuf_r);CHKERRQ(ierr);
13720572522cSBarry Smith   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
137338571addSHong Zhang #if defined(PTAP_PROFILE)
13748563dfccSBarry Smith   ierr = PetscTime(&t4);CHKERRQ(ierr);
13759516a85cSHong Zhang   if (rank==1) {
137605d62848SHong 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);
13779516a85cSHong Zhang   }
137838571addSHong Zhang #endif
137942c57489SHong Zhang   PetscFunctionReturn(0);
138042c57489SHong Zhang }
1381