xref: /petsc/src/mat/impls/aij/mpi/mpimatmatmult.c (revision 33ba5e2fdab343b05366fddedfd8875376f7a774)
1be1d678aSKris Buschelman 
22499ec78SHong Zhang /*
32499ec78SHong Zhang   Defines matrix-matrix product routines for pairs of MPIAIJ matrices
42499ec78SHong Zhang           C = A * B
52499ec78SHong Zhang */
6c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/aij.h> /*I "petscmat.h" I*/
7c6db04a5SJed Brown #include <../src/mat/utils/freespace.h>
8c6db04a5SJed Brown #include <../src/mat/impls/aij/mpi/mpiaij.h>
9c6db04a5SJed Brown #include <petscbt.h>
10c6db04a5SJed Brown #include <../src/mat/impls/dense/mpi/mpidense.h>
11af0996ceSBarry Smith #include <petsc/private/vecimpl.h>
121634b032SHong Zhang 
132499ec78SHong Zhang #undef __FUNCT__
142499ec78SHong Zhang #define __FUNCT__ "MatMatMult_MPIAIJ_MPIAIJ"
15150d2497SBarry Smith PETSC_INTERN PetscErrorCode MatMatMult_MPIAIJ_MPIAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill, Mat *C)
162499ec78SHong Zhang {
172499ec78SHong Zhang   PetscErrorCode ierr;
180fc8cf34SHong Zhang   const char     *algTypes[2] = {"scalable","nonscalable"};
190d3441aeSHong Zhang   PetscInt       alg=1; /* set default algorithm */
20e55be12dSHong Zhang   MPI_Comm       comm;
212499ec78SHong Zhang 
222499ec78SHong Zhang   PetscFunctionBegin;
232499ec78SHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
24e55be12dSHong Zhang     ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
256c4ed002SBarry Smith     if (A->cmap->rstart != B->rmap->rstart || A->cmap->rend != B->rmap->rend) SETERRQ4(comm,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, (%D, %D) != (%D,%D)",A->cmap->rstart,A->cmap->rend,B->rmap->rstart,B->rmap->rend);
26e55be12dSHong Zhang 
270fc8cf34SHong Zhang     ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr);
2878d30b94SHong Zhang     ierr = PetscOptionsEList("-matmatmult_via","Algorithmic approach","MatMatMult",algTypes,2,algTypes[1],&alg,NULL);CHKERRQ(ierr);
290fc8cf34SHong Zhang     ierr = PetscOptionsEnd();CHKERRQ(ierr);
300fc8cf34SHong Zhang 
313ff4c91cSHong Zhang     ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
320fc8cf34SHong Zhang     switch (alg) {
330fc8cf34SHong Zhang     case 1:
340fc8cf34SHong Zhang       ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(A,B,fill,C);CHKERRQ(ierr);
350fc8cf34SHong Zhang       break;
360fc8cf34SHong Zhang     default:
37b2405163SHong Zhang       ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ(A,B,fill,C);CHKERRQ(ierr);
380fc8cf34SHong Zhang       break;
390fc8cf34SHong Zhang     }
403ff4c91cSHong Zhang     ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
412499ec78SHong Zhang   }
423ff4c91cSHong Zhang   ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
43598bc09dSHong Zhang   ierr = (*(*C)->ops->matmultnumeric)(A,B,*C);CHKERRQ(ierr);
443ff4c91cSHong Zhang   ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
452499ec78SHong Zhang   PetscFunctionReturn(0);
462499ec78SHong Zhang }
472499ec78SHong Zhang 
48f33d1a9aSHong Zhang #undef __FUNCT__
49a1a4e84aSHong Zhang #define __FUNCT__ "MatDestroy_MPIAIJ_MatMatMult"
50a1a4e84aSHong Zhang PetscErrorCode MatDestroy_MPIAIJ_MatMatMult(Mat A)
51598bc09dSHong Zhang {
52598bc09dSHong Zhang   PetscErrorCode ierr;
53598bc09dSHong Zhang   Mat_MPIAIJ     *a    = (Mat_MPIAIJ*)A->data;
54598bc09dSHong Zhang   Mat_PtAPMPI    *ptap = a->ptap;
55598bc09dSHong Zhang 
56598bc09dSHong Zhang   PetscFunctionBegin;
57b7f45c76SHong Zhang   ierr = PetscFree2(ptap->startsj_s,ptap->startsj_r);CHKERRQ(ierr);
58598bc09dSHong Zhang   ierr = PetscFree(ptap->bufa);CHKERRQ(ierr);
59a1a4e84aSHong Zhang   ierr = MatDestroy(&ptap->P_loc);CHKERRQ(ierr);
60a1a4e84aSHong Zhang   ierr = MatDestroy(&ptap->P_oth);CHKERRQ(ierr);
61ab1b013aSHong Zhang   ierr = MatDestroy(&ptap->Pt);CHKERRQ(ierr);
62a1a4e84aSHong Zhang   ierr = PetscFree(ptap->api);CHKERRQ(ierr);
63a1a4e84aSHong Zhang   ierr = PetscFree(ptap->apj);CHKERRQ(ierr);
64598bc09dSHong Zhang   ierr = PetscFree(ptap->apa);CHKERRQ(ierr);
65598bc09dSHong Zhang   ierr = ptap->destroy(A);CHKERRQ(ierr);
66598bc09dSHong Zhang   ierr = PetscFree(ptap);CHKERRQ(ierr);
67598bc09dSHong Zhang   PetscFunctionReturn(0);
68598bc09dSHong Zhang }
69598bc09dSHong Zhang 
702499ec78SHong Zhang #undef __FUNCT__
71a1a4e84aSHong Zhang #define __FUNCT__ "MatDuplicate_MPIAIJ_MatMatMult"
72a1a4e84aSHong Zhang PetscErrorCode MatDuplicate_MPIAIJ_MatMatMult(Mat A, MatDuplicateOption op, Mat *M)
734ae0847bSHong Zhang {
744ae0847bSHong Zhang   PetscErrorCode ierr;
754ae0847bSHong Zhang   Mat_MPIAIJ     *a    = (Mat_MPIAIJ*)A->data;
764ae0847bSHong Zhang   Mat_PtAPMPI    *ptap = a->ptap;
774ae0847bSHong Zhang 
784ae0847bSHong Zhang   PetscFunctionBegin;
794ae0847bSHong Zhang   ierr = (*ptap->duplicate)(A,op,M);CHKERRQ(ierr);
8026fbe8dcSKarl Rupp 
814ae0847bSHong Zhang   (*M)->ops->destroy   = ptap->destroy;   /* = MatDestroy_MPIAIJ, *M doesn't duplicate A's special structure! */
824ae0847bSHong Zhang   (*M)->ops->duplicate = ptap->duplicate; /* = MatDuplicate_MPIAIJ */
834ae0847bSHong Zhang   PetscFunctionReturn(0);
844ae0847bSHong Zhang }
854ae0847bSHong Zhang 
864ae0847bSHong Zhang #undef __FUNCT__
87f8efd71cSHong Zhang #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable"
88f8efd71cSHong Zhang PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable(Mat A,Mat P,Mat C)
89598bc09dSHong Zhang {
90598bc09dSHong Zhang   PetscErrorCode ierr;
914ae0847bSHong Zhang   Mat_MPIAIJ     *a  =(Mat_MPIAIJ*)A->data,*c=(Mat_MPIAIJ*)C->data;
92598bc09dSHong Zhang   Mat_SeqAIJ     *ad =(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
93598bc09dSHong Zhang   Mat_SeqAIJ     *cd =(Mat_SeqAIJ*)(c->A)->data,*co=(Mat_SeqAIJ*)(c->B)->data;
94c12557a1SHong Zhang   PetscScalar    *cda=cd->a,*coa=co->a;
95598bc09dSHong Zhang   Mat_SeqAIJ     *p_loc,*p_oth;
969ce11a7cSHong Zhang   PetscScalar    *apa,*ca;
97c12557a1SHong Zhang   PetscInt       cm   =C->rmap->n;
98598bc09dSHong Zhang   Mat_PtAPMPI    *ptap=c->ptap;
99c12557a1SHong Zhang   PetscInt       *api,*apj,*apJ,i,k;
10029825010SHong Zhang   PetscInt       cstart=C->cmap->rstart;
101598bc09dSHong Zhang   PetscInt       cdnz,conz,k0,k1;
1029467f45fSHong Zhang   MPI_Comm       comm;
1039467f45fSHong Zhang   PetscMPIInt    size;
104598bc09dSHong Zhang 
105598bc09dSHong Zhang   PetscFunctionBegin;
1069467f45fSHong Zhang   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
1079467f45fSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1089467f45fSHong Zhang 
109a1a4e84aSHong Zhang   /* 1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
110598bc09dSHong Zhang   /*-----------------------------------------------------*/
111a1a4e84aSHong Zhang   /* update numerical values of P_oth and P_loc */
112b7f45c76SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
113a1a4e84aSHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
114598bc09dSHong Zhang 
115598bc09dSHong Zhang   /* 2) compute numeric C_loc = A_loc*P = Ad*P_loc + Ao*P_oth */
116598bc09dSHong Zhang   /*----------------------------------------------------------*/
117598bc09dSHong Zhang   /* get data from symbolic products */
118a1a4e84aSHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
119c12557a1SHong Zhang   p_oth = NULL;
1209467f45fSHong Zhang   if (size >1) {
1219467f45fSHong Zhang     p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
1229467f45fSHong Zhang   }
123598bc09dSHong Zhang 
124598bc09dSHong Zhang   /* get apa for storing dense row A[i,:]*P */
125598bc09dSHong Zhang   apa = ptap->apa;
126598bc09dSHong Zhang 
12729825010SHong Zhang   api = ptap->api;
12829825010SHong Zhang   apj = ptap->apj;
129598bc09dSHong Zhang   for (i=0; i<cm; i++) {
130c12557a1SHong Zhang     /* compute apa = A[i,:]*P */
131e497d3c8SHong Zhang     AProw_nonscalable(i,ad,ao,p_loc,p_oth,apa);
132598bc09dSHong Zhang 
133598bc09dSHong Zhang     /* set values in C */
134598bc09dSHong Zhang     apJ  = apj + api[i];
135598bc09dSHong Zhang     cdnz = cd->i[i+1] - cd->i[i];
136598bc09dSHong Zhang     conz = co->i[i+1] - co->i[i];
137598bc09dSHong Zhang 
138598bc09dSHong Zhang     /* 1st off-diagoanl part of C */
139598bc09dSHong Zhang     ca = coa + co->i[i];
140598bc09dSHong Zhang     k  = 0;
141598bc09dSHong Zhang     for (k0=0; k0<conz; k0++) {
142598bc09dSHong Zhang       if (apJ[k] >= cstart) break;
143598bc09dSHong Zhang       ca[k0]      = apa[apJ[k]];
144598bc09dSHong Zhang       apa[apJ[k]] = 0.0;
145598bc09dSHong Zhang       k++;
146598bc09dSHong Zhang     }
147598bc09dSHong Zhang 
148598bc09dSHong Zhang     /* diagonal part of C */
149598bc09dSHong Zhang     ca = cda + cd->i[i];
150598bc09dSHong Zhang     for (k1=0; k1<cdnz; k1++) {
151598bc09dSHong Zhang       ca[k1]      = apa[apJ[k]];
152598bc09dSHong Zhang       apa[apJ[k]] = 0.0;
153598bc09dSHong Zhang       k++;
154598bc09dSHong Zhang     }
155598bc09dSHong Zhang 
156598bc09dSHong Zhang     /* 2nd off-diagoanl part of C */
157598bc09dSHong Zhang     ca = coa + co->i[i];
158598bc09dSHong Zhang     for (; k0<conz; k0++) {
159598bc09dSHong Zhang       ca[k0]      = apa[apJ[k]];
160598bc09dSHong Zhang       apa[apJ[k]] = 0.0;
161598bc09dSHong Zhang       k++;
162598bc09dSHong Zhang     }
163598bc09dSHong Zhang   }
164598bc09dSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
165598bc09dSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
166598bc09dSHong Zhang   PetscFunctionReturn(0);
167598bc09dSHong Zhang }
168598bc09dSHong Zhang 
169598bc09dSHong Zhang #undef __FUNCT__
1700fc8cf34SHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable"
1710fc8cf34SHong Zhang PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat A,Mat P,PetscReal fill,Mat *C)
172598bc09dSHong Zhang {
173598bc09dSHong Zhang   PetscErrorCode     ierr;
174ce94432eSBarry Smith   MPI_Comm           comm;
1759467f45fSHong Zhang   PetscMPIInt        size;
176bfcd1a73SHong Zhang   Mat                Cmpi;
177598bc09dSHong Zhang   Mat_PtAPMPI        *ptap;
1780298fd71SBarry Smith   PetscFreeSpaceList free_space=NULL,current_space=NULL;
1794ae0847bSHong Zhang   Mat_MPIAIJ         *a        =(Mat_MPIAIJ*)A->data,*c;
180bfcd1a73SHong Zhang   Mat_SeqAIJ         *ad       =(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_loc,*p_oth;
1814ae0847bSHong Zhang   PetscInt           *pi_loc,*pj_loc,*pi_oth,*pj_oth,*dnz,*onz;
1824ae0847bSHong Zhang   PetscInt           *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,rstart=A->rmap->rstart;
183d14fa243SHong Zhang   PetscInt           *lnk,i,pnz,row,*api,*apj,*Jptr,apnz,nspacedouble=0,j,nzi;
184e55be12dSHong Zhang   PetscInt           am=A->rmap->n,pN=P->cmap->N,pn=P->cmap->n,pm=P->rmap->n,Crmax;
185598bc09dSHong Zhang   PetscBT            lnkbt;
186598bc09dSHong Zhang   PetscScalar        *apa;
187bfcd1a73SHong Zhang   PetscReal          afill;
188e55be12dSHong Zhang   PetscTable         ta;
189598bc09dSHong Zhang 
190598bc09dSHong Zhang   PetscFunctionBegin;
191ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
1929467f45fSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1939467f45fSHong Zhang 
194598bc09dSHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
195b00a9115SJed Brown   ierr = PetscNew(&ptap);CHKERRQ(ierr);
196598bc09dSHong Zhang 
197598bc09dSHong Zhang   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
198b7f45c76SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
19919f0e57aSHong Zhang 
200598bc09dSHong Zhang   /* get P_loc by taking all local rows of P */
201a1a4e84aSHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
202598bc09dSHong Zhang 
203a1a4e84aSHong Zhang   p_loc  = (Mat_SeqAIJ*)(ptap->P_loc)->data;
204598bc09dSHong Zhang   pi_loc = p_loc->i; pj_loc = p_loc->j;
2059467f45fSHong Zhang   if (size > 1) {
2069467f45fSHong Zhang     p_oth  = (Mat_SeqAIJ*)(ptap->P_oth)->data;
207598bc09dSHong Zhang     pi_oth = p_oth->i; pj_oth = p_oth->j;
20820e3dc75SHong Zhang   } else {
20920e3dc75SHong Zhang     p_oth = NULL;
21020e3dc75SHong Zhang     pi_oth = NULL; pj_oth = NULL;
2119467f45fSHong Zhang   }
212598bc09dSHong Zhang 
213598bc09dSHong Zhang   /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */
214598bc09dSHong Zhang   /*-------------------------------------------------------------------*/
215854ce69bSBarry Smith   ierr      = PetscMalloc1(am+2,&api);CHKERRQ(ierr);
216a1a4e84aSHong Zhang   ptap->api = api;
217598bc09dSHong Zhang   api[0]    = 0;
218598bc09dSHong Zhang 
21945d00d1dSHong Zhang   /* create and initialize a linked list -- TODO: replace it with PetscBTCreate()! */
22045d00d1dSHong Zhang   ierr = PetscTableCreate(pn,pN,&ta);CHKERRQ(ierr);
2214b38b95cSHong Zhang   MatRowMergeMax_SeqAIJ(p_loc,ptap->P_loc->rmap->N,ta);
2224b38b95cSHong Zhang   MatRowMergeMax_SeqAIJ(p_oth,ptap->P_oth->rmap->N,ta);
223e55be12dSHong Zhang   ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr);
224e55be12dSHong Zhang   ierr = PetscTableDestroy(&ta);CHKERRQ(ierr);
225e55be12dSHong Zhang 
226e55be12dSHong Zhang   ierr = PetscLLCondensedCreate(Crmax,pN,&lnk,&lnkbt);CHKERRQ(ierr);
227598bc09dSHong Zhang 
228bfcd1a73SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */
229f91af8c7SBarry Smith   ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(adi[am],PetscIntSumTruncate(aoi[am],pi_loc[pm]))),&free_space);CHKERRQ(ierr);
230598bc09dSHong Zhang   current_space = free_space;
231598bc09dSHong Zhang 
232598bc09dSHong Zhang   ierr = MatPreallocateInitialize(comm,am,pn,dnz,onz);CHKERRQ(ierr);
233598bc09dSHong Zhang   for (i=0; i<am; i++) {
234598bc09dSHong Zhang     /* diagonal portion of A */
235598bc09dSHong Zhang     nzi = adi[i+1] - adi[i];
236598bc09dSHong Zhang     for (j=0; j<nzi; j++) {
237598bc09dSHong Zhang       row  = *adj++;
238598bc09dSHong Zhang       pnz  = pi_loc[row+1] - pi_loc[row];
239598bc09dSHong Zhang       Jptr = pj_loc + pi_loc[row];
240598bc09dSHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
2411fbbb359SHong Zhang       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
242598bc09dSHong Zhang     }
243598bc09dSHong Zhang     /* off-diagonal portion of A */
244598bc09dSHong Zhang     nzi = aoi[i+1] - aoi[i];
245598bc09dSHong Zhang     for (j=0; j<nzi; j++) {
246598bc09dSHong Zhang       row  = *aoj++;
247598bc09dSHong Zhang       pnz  = pi_oth[row+1] - pi_oth[row];
248598bc09dSHong Zhang       Jptr = pj_oth + pi_oth[row];
2491fbbb359SHong Zhang       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
250598bc09dSHong Zhang     }
251598bc09dSHong Zhang 
252d14fa243SHong Zhang     apnz     = lnk[0];
253598bc09dSHong Zhang     api[i+1] = api[i] + apnz;
254598bc09dSHong Zhang 
255598bc09dSHong Zhang     /* if free space is not available, double the total space in the list */
256598bc09dSHong Zhang     if (current_space->local_remaining<apnz) {
257f91af8c7SBarry Smith       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
258598bc09dSHong Zhang       nspacedouble++;
259598bc09dSHong Zhang     }
260598bc09dSHong Zhang 
261598bc09dSHong Zhang     /* Copy data into free space, then initialize lnk */
2621fbbb359SHong Zhang     ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
263598bc09dSHong Zhang     ierr = MatPreallocateSet(i+rstart,apnz,current_space->array,dnz,onz);CHKERRQ(ierr);
2642205254eSKarl Rupp 
265598bc09dSHong Zhang     current_space->array           += apnz;
266598bc09dSHong Zhang     current_space->local_used      += apnz;
267598bc09dSHong Zhang     current_space->local_remaining -= apnz;
268598bc09dSHong Zhang   }
269598bc09dSHong Zhang 
270598bc09dSHong Zhang   /* Allocate space for apj, initialize apj, and */
271598bc09dSHong Zhang   /* destroy list of free space and other temporary array(s) */
272854ce69bSBarry Smith   ierr = PetscMalloc1(api[am]+1,&ptap->apj);CHKERRQ(ierr);
273a1a4e84aSHong Zhang   apj  = ptap->apj;
274a1a4e84aSHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,ptap->apj);CHKERRQ(ierr);
275598bc09dSHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
276598bc09dSHong Zhang 
277f84c536bSHong Zhang   /* malloc apa to store dense row A[i,:]*P */
2781795a4d1SJed Brown   ierr = PetscCalloc1(pN,&apa);CHKERRQ(ierr);
2792205254eSKarl Rupp 
280f84c536bSHong Zhang   ptap->apa = apa;
281f84c536bSHong Zhang 
282bfcd1a73SHong Zhang   /* create and assemble symbolic parallel matrix Cmpi */
283598bc09dSHong Zhang   /*----------------------------------------------------*/
284bfcd1a73SHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
285bfcd1a73SHong Zhang   ierr = MatSetSizes(Cmpi,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
28633d57670SJed Brown   ierr = MatSetBlockSizesFromMats(Cmpi,A,P);CHKERRQ(ierr);
287a2f3521dSMark F. Adams 
288bfcd1a73SHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
289bfcd1a73SHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
290598bc09dSHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
291598bc09dSHong Zhang   for (i=0; i<am; i++) {
292598bc09dSHong Zhang     row  = i + rstart;
293598bc09dSHong Zhang     apnz = api[i+1] - api[i];
294bfcd1a73SHong Zhang     ierr = MatSetValues(Cmpi,1,&row,apnz,apj,apa,INSERT_VALUES);CHKERRQ(ierr);
295598bc09dSHong Zhang     apj += apnz;
296598bc09dSHong Zhang   }
297bfcd1a73SHong Zhang   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
298bfcd1a73SHong Zhang   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
299598bc09dSHong Zhang 
300bfcd1a73SHong Zhang   ptap->destroy        = Cmpi->ops->destroy;
301bfcd1a73SHong Zhang   ptap->duplicate      = Cmpi->ops->duplicate;
302f8efd71cSHong Zhang   Cmpi->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable;
303bfcd1a73SHong Zhang   Cmpi->ops->destroy   = MatDestroy_MPIAIJ_MatMatMult;
304bfcd1a73SHong Zhang   Cmpi->ops->duplicate = MatDuplicate_MPIAIJ_MatMatMult;
305598bc09dSHong Zhang 
306bfcd1a73SHong Zhang   /* attach the supporting struct to Cmpi for reuse */
307bfcd1a73SHong Zhang   c       = (Mat_MPIAIJ*)Cmpi->data;
308598bc09dSHong Zhang   c->ptap = ptap;
309598bc09dSHong Zhang 
310bfcd1a73SHong Zhang   *C = Cmpi;
311bfcd1a73SHong Zhang 
312bfcd1a73SHong Zhang   /* set MatInfo */
313118727c9SMark F. Adams   afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1) + 1.e-5;
314bfcd1a73SHong Zhang   if (afill < 1.0) afill = 1.0;
315bfcd1a73SHong Zhang   Cmpi->info.mallocs           = nspacedouble;
316bfcd1a73SHong Zhang   Cmpi->info.fill_ratio_given  = fill;
317bfcd1a73SHong Zhang   Cmpi->info.fill_ratio_needed = afill;
318bfcd1a73SHong Zhang 
319bfcd1a73SHong Zhang #if defined(PETSC_USE_INFO)
320bfcd1a73SHong Zhang   if (api[am]) {
32157622a8eSBarry Smith     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
32257622a8eSBarry Smith     ierr = PetscInfo1(Cmpi,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr);
323bfcd1a73SHong Zhang   } else {
324bfcd1a73SHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
325bfcd1a73SHong Zhang   }
326bfcd1a73SHong Zhang #endif
327598bc09dSHong Zhang   PetscFunctionReturn(0);
328598bc09dSHong Zhang }
329598bc09dSHong Zhang 
3309123193aSHong Zhang #undef __FUNCT__
3319123193aSHong Zhang #define __FUNCT__ "MatMatMult_MPIAIJ_MPIDense"
332150d2497SBarry Smith PETSC_INTERN PetscErrorCode MatMatMult_MPIAIJ_MPIDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
3339123193aSHong Zhang {
3349123193aSHong Zhang   PetscErrorCode ierr;
3359123193aSHong Zhang 
3369123193aSHong Zhang   PetscFunctionBegin;
3379123193aSHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
3383ff4c91cSHong Zhang     ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
3399123193aSHong Zhang     ierr = MatMatMultSymbolic_MPIAIJ_MPIDense(A,B,fill,C);CHKERRQ(ierr);
3403ff4c91cSHong Zhang     ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
3419123193aSHong Zhang   }
3423ff4c91cSHong Zhang   ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
3439123193aSHong Zhang   ierr = MatMatMultNumeric_MPIAIJ_MPIDense(A,B,*C);CHKERRQ(ierr);
3443ff4c91cSHong Zhang   ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
3459123193aSHong Zhang   PetscFunctionReturn(0);
3469123193aSHong Zhang }
3479123193aSHong Zhang 
3484324174eSBarry Smith typedef struct {
3494324174eSBarry Smith   Mat         workB;
3504324174eSBarry Smith   PetscScalar *rvalues,*svalues;
3514324174eSBarry Smith   MPI_Request *rwaits,*swaits;
3524324174eSBarry Smith } MPIAIJ_MPIDense;
3534324174eSBarry Smith 
3547af9e4fdSHong Zhang #undef __FUNCT__
35596b95a6bSBarry Smith #define __FUNCT__ "MatMPIAIJ_MPIDenseDestroy"
35696b95a6bSBarry Smith PetscErrorCode MatMPIAIJ_MPIDenseDestroy(void *ctx)
3574324174eSBarry Smith {
3584324174eSBarry Smith   MPIAIJ_MPIDense *contents = (MPIAIJ_MPIDense*) ctx;
3594324174eSBarry Smith   PetscErrorCode  ierr;
3604324174eSBarry Smith 
3614324174eSBarry Smith   PetscFunctionBegin;
3626bf464f9SBarry Smith   ierr = MatDestroy(&contents->workB);CHKERRQ(ierr);
3634324174eSBarry Smith   ierr = PetscFree4(contents->rvalues,contents->svalues,contents->rwaits,contents->swaits);CHKERRQ(ierr);
3644324174eSBarry Smith   ierr = PetscFree(contents);CHKERRQ(ierr);
3654324174eSBarry Smith   PetscFunctionReturn(0);
3664324174eSBarry Smith }
3674324174eSBarry Smith 
3689123193aSHong Zhang #undef __FUNCT__
369fd4e9aacSBarry Smith #define __FUNCT__ "MatMatMultNumeric_MPIDense"
370fd4e9aacSBarry Smith /*
371fd4e9aacSBarry Smith     This is a "dummy function" that handles the case where matrix C was created as a dense matrix
372fd4e9aacSBarry Smith   directly by the user and passed to MatMatMult() with the MAT_REUSE_MATRIX option
373fd4e9aacSBarry Smith 
374fd4e9aacSBarry Smith   It is the same as MatMatMultSymbolic_MPIAIJ_MPIDense() except does not create C
375fd4e9aacSBarry Smith */
376fd4e9aacSBarry Smith PetscErrorCode MatMatMultNumeric_MPIDense(Mat A,Mat B,Mat C)
377fd4e9aacSBarry Smith {
378fd4e9aacSBarry Smith   PetscErrorCode         ierr;
379fd4e9aacSBarry Smith   PetscBool              flg;
380fd4e9aacSBarry Smith   Mat_MPIAIJ             *aij = (Mat_MPIAIJ*) A->data;
381fd4e9aacSBarry Smith   PetscInt               nz   = aij->B->cmap->n;
382fd4e9aacSBarry Smith   PetscContainer         container;
383fd4e9aacSBarry Smith   MPIAIJ_MPIDense        *contents;
384fd4e9aacSBarry Smith   VecScatter             ctx   = aij->Mvctx;
385fd4e9aacSBarry Smith   VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata;
386fd4e9aacSBarry Smith   VecScatter_MPI_General *to   = (VecScatter_MPI_General*) ctx->todata;
387fd4e9aacSBarry Smith 
388fd4e9aacSBarry Smith   PetscFunctionBegin;
389fd4e9aacSBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)B,MATMPIDENSE,&flg);CHKERRQ(ierr);
390fd4e9aacSBarry Smith   if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Second matrix must be mpidense");
391fd4e9aacSBarry Smith 
392fd4e9aacSBarry Smith   /* Handle case where where user provided the final C matrix rather than calling MatMatMult() with MAT_INITIAL_MATRIX*/
393fd4e9aacSBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&flg);CHKERRQ(ierr);
394fd4e9aacSBarry Smith   if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"First matrix must be MPIAIJ");
395fd4e9aacSBarry Smith 
396fd4e9aacSBarry Smith   C->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIDense;
397fd4e9aacSBarry Smith 
398fd4e9aacSBarry Smith   ierr = PetscNew(&contents);CHKERRQ(ierr);
399fd4e9aacSBarry Smith   /* Create work matrix used to store off processor rows of B needed for local product */
400fd4e9aacSBarry Smith   ierr = MatCreateSeqDense(PETSC_COMM_SELF,nz,B->cmap->N,NULL,&contents->workB);CHKERRQ(ierr);
401fd4e9aacSBarry Smith   /* Create work arrays needed */
402fd4e9aacSBarry Smith   ierr = PetscMalloc4(B->cmap->N*from->starts[from->n],&contents->rvalues,
403fd4e9aacSBarry Smith                       B->cmap->N*to->starts[to->n],&contents->svalues,
404fd4e9aacSBarry Smith                       from->n,&contents->rwaits,
405fd4e9aacSBarry Smith                       to->n,&contents->swaits);CHKERRQ(ierr);
406fd4e9aacSBarry Smith 
407fd4e9aacSBarry Smith   ierr = PetscContainerCreate(PetscObjectComm((PetscObject)A),&container);CHKERRQ(ierr);
408fd4e9aacSBarry Smith   ierr = PetscContainerSetPointer(container,contents);CHKERRQ(ierr);
409fd4e9aacSBarry Smith   ierr = PetscContainerSetUserDestroy(container,MatMPIAIJ_MPIDenseDestroy);CHKERRQ(ierr);
410fd4e9aacSBarry Smith   ierr = PetscObjectCompose((PetscObject)C,"workB",(PetscObject)container);CHKERRQ(ierr);
411fd4e9aacSBarry Smith   ierr = PetscContainerDestroy(&container);CHKERRQ(ierr);
412fd4e9aacSBarry Smith 
413fd4e9aacSBarry Smith   ierr = (*C->ops->matmultnumeric)(A,B,C);CHKERRQ(ierr);
414fd4e9aacSBarry Smith   PetscFunctionReturn(0);
415fd4e9aacSBarry Smith }
416fd4e9aacSBarry Smith 
417fd4e9aacSBarry Smith #undef __FUNCT__
4189123193aSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIDense"
4199123193aSHong Zhang PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIDense(Mat A,Mat B,PetscReal fill,Mat *C)
4209123193aSHong Zhang {
4219123193aSHong Zhang   PetscErrorCode         ierr;
422f32d5d43SBarry Smith   Mat_MPIAIJ             *aij = (Mat_MPIAIJ*) A->data;
423d0f46423SBarry Smith   PetscInt               nz   = aij->B->cmap->n;
424bf0cc555SLisandro Dalcin   PetscContainer         container;
4254324174eSBarry Smith   MPIAIJ_MPIDense        *contents;
4264324174eSBarry Smith   VecScatter             ctx   = aij->Mvctx;
4274324174eSBarry Smith   VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata;
4284324174eSBarry Smith   VecScatter_MPI_General *to   = (VecScatter_MPI_General*) ctx->todata;
429d0f46423SBarry Smith   PetscInt               m     = A->rmap->n,n=B->cmap->n;
4309123193aSHong Zhang 
4319123193aSHong Zhang   PetscFunctionBegin;
432ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)B),C);CHKERRQ(ierr);
433d0f46423SBarry Smith   ierr = MatSetSizes(*C,m,n,A->rmap->N,B->cmap->N);CHKERRQ(ierr);
43433d57670SJed Brown   ierr = MatSetBlockSizesFromMats(*C,A,B);CHKERRQ(ierr);
435cb2480eaSBarry Smith   ierr = MatSetType(*C,MATMPIDENSE);CHKERRQ(ierr);
4360298fd71SBarry Smith   ierr = MatMPIDenseSetPreallocation(*C,NULL);CHKERRQ(ierr);
437cb2480eaSBarry Smith   ierr = MatAssemblyBegin(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
438cb2480eaSBarry Smith   ierr = MatAssemblyEnd(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
4392205254eSKarl Rupp 
440d73949e8SHong Zhang   (*C)->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIDense;
441f32d5d43SBarry Smith 
442b00a9115SJed Brown   ierr = PetscNew(&contents);CHKERRQ(ierr);
443f32d5d43SBarry Smith   /* Create work matrix used to store off processor rows of B needed for local product */
4440298fd71SBarry Smith   ierr = MatCreateSeqDense(PETSC_COMM_SELF,nz,B->cmap->N,NULL,&contents->workB);CHKERRQ(ierr);
4454324174eSBarry Smith   /* Create work arrays needed */
446dcca6d9dSJed Brown   ierr = PetscMalloc4(B->cmap->N*from->starts[from->n],&contents->rvalues,
447dcca6d9dSJed Brown                       B->cmap->N*to->starts[to->n],&contents->svalues,
448dcca6d9dSJed Brown                       from->n,&contents->rwaits,
449dcca6d9dSJed Brown                       to->n,&contents->swaits);CHKERRQ(ierr);
4504324174eSBarry Smith 
451ce94432eSBarry Smith   ierr = PetscContainerCreate(PetscObjectComm((PetscObject)A),&container);CHKERRQ(ierr);
452bf0cc555SLisandro Dalcin   ierr = PetscContainerSetPointer(container,contents);CHKERRQ(ierr);
45396b95a6bSBarry Smith   ierr = PetscContainerSetUserDestroy(container,MatMPIAIJ_MPIDenseDestroy);CHKERRQ(ierr);
454bf0cc555SLisandro Dalcin   ierr = PetscObjectCompose((PetscObject)(*C),"workB",(PetscObject)container);CHKERRQ(ierr);
455bf0cc555SLisandro Dalcin   ierr = PetscContainerDestroy(&container);CHKERRQ(ierr);
4569123193aSHong Zhang   PetscFunctionReturn(0);
4579123193aSHong Zhang }
4589123193aSHong Zhang 
4597af9e4fdSHong Zhang #undef __FUNCT__
4607af9e4fdSHong Zhang #define __FUNCT__ "MatMPIDenseScatter"
461f32d5d43SBarry Smith /*
4622636ff26SBarry Smith     Performs an efficient scatter on the rows of B needed by this process; this is
4632636ff26SBarry Smith     a modification of the VecScatterBegin_() routines.
464f32d5d43SBarry Smith */
4654324174eSBarry Smith PetscErrorCode MatMPIDenseScatter(Mat A,Mat B,Mat C,Mat *outworkB)
466f32d5d43SBarry Smith {
467f32d5d43SBarry Smith   Mat_MPIAIJ             *aij = (Mat_MPIAIJ*)A->data;
468f32d5d43SBarry Smith   PetscErrorCode         ierr;
469f32d5d43SBarry Smith   PetscScalar            *b,*w,*svalues,*rvalues;
470f32d5d43SBarry Smith   VecScatter             ctx   = aij->Mvctx;
471f32d5d43SBarry Smith   VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata;
472f32d5d43SBarry Smith   VecScatter_MPI_General *to   = (VecScatter_MPI_General*) ctx->todata;
473f32d5d43SBarry Smith   PetscInt               i,j,k;
474aa5bb8c0SSatish Balay   PetscInt               *sindices,*sstarts,*rindices,*rstarts;
475aa5bb8c0SSatish Balay   PetscMPIInt            *sprocs,*rprocs,nrecvs;
476f32d5d43SBarry Smith   MPI_Request            *swaits,*rwaits;
477ce94432eSBarry Smith   MPI_Comm               comm;
478d0f46423SBarry Smith   PetscMPIInt            tag  = ((PetscObject)ctx)->tag,ncols = B->cmap->N, nrows = aij->B->cmap->n,imdex,nrowsB = B->rmap->n;
479f32d5d43SBarry Smith   MPI_Status             status;
4804324174eSBarry Smith   MPIAIJ_MPIDense        *contents;
481bf0cc555SLisandro Dalcin   PetscContainer         container;
4824324174eSBarry Smith   Mat                    workB;
483f32d5d43SBarry Smith 
484f32d5d43SBarry Smith   PetscFunctionBegin;
485ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
486bf0cc555SLisandro Dalcin   ierr = PetscObjectQuery((PetscObject)C,"workB",(PetscObject*)&container);CHKERRQ(ierr);
48729825010SHong Zhang   if (!container) SETERRQ(comm,PETSC_ERR_PLIB,"Container does not exist");
488bf0cc555SLisandro Dalcin   ierr = PetscContainerGetPointer(container,(void**)&contents);CHKERRQ(ierr);
4894324174eSBarry Smith 
4904324174eSBarry Smith   workB = *outworkB = contents->workB;
491cf1a0672SHong Zhang   if (nrows != workB->rmap->n) SETERRQ2(comm,PETSC_ERR_PLIB,"Number of rows of workB %D not equal to columns of aij->B %D",nrows,workB->cmap->n);
492f32d5d43SBarry Smith   sindices = to->indices;
493f32d5d43SBarry Smith   sstarts  = to->starts;
494f32d5d43SBarry Smith   sprocs   = to->procs;
4954324174eSBarry Smith   swaits   = contents->swaits;
4964324174eSBarry Smith   svalues  = contents->svalues;
497f32d5d43SBarry Smith 
498f32d5d43SBarry Smith   rindices = from->indices;
499f32d5d43SBarry Smith   rstarts  = from->starts;
500f32d5d43SBarry Smith   rprocs   = from->procs;
5014324174eSBarry Smith   rwaits   = contents->rwaits;
5024324174eSBarry Smith   rvalues  = contents->rvalues;
503f32d5d43SBarry Smith 
5048c778c55SBarry Smith   ierr = MatDenseGetArray(B,&b);CHKERRQ(ierr);
5058c778c55SBarry Smith   ierr = MatDenseGetArray(workB,&w);CHKERRQ(ierr);
506f32d5d43SBarry Smith 
507f32d5d43SBarry Smith   for (i=0; i<from->n; i++) {
508f32d5d43SBarry Smith     ierr = MPI_Irecv(rvalues+ncols*rstarts[i],ncols*(rstarts[i+1]-rstarts[i]),MPIU_SCALAR,rprocs[i],tag,comm,rwaits+i);CHKERRQ(ierr);
509f32d5d43SBarry Smith   }
510f32d5d43SBarry Smith 
511f32d5d43SBarry Smith   for (i=0; i<to->n; i++) {
512f32d5d43SBarry Smith     /* pack a message at a time */
513f32d5d43SBarry Smith     for (j=0; j<sstarts[i+1]-sstarts[i]; j++) {
514f32d5d43SBarry Smith       for (k=0; k<ncols; k++) {
515f32d5d43SBarry Smith         svalues[ncols*(sstarts[i] + j) + k] = b[sindices[sstarts[i]+j] + nrowsB*k];
5162636ff26SBarry Smith       }
5172636ff26SBarry Smith     }
518f32d5d43SBarry Smith     ierr = MPI_Isend(svalues+ncols*sstarts[i],ncols*(sstarts[i+1]-sstarts[i]),MPIU_SCALAR,sprocs[i],tag,comm,swaits+i);CHKERRQ(ierr);
519f32d5d43SBarry Smith   }
5202636ff26SBarry Smith 
521f32d5d43SBarry Smith   nrecvs = from->n;
522f32d5d43SBarry Smith   while (nrecvs) {
523f32d5d43SBarry Smith     ierr = MPI_Waitany(from->n,rwaits,&imdex,&status);CHKERRQ(ierr);
524f32d5d43SBarry Smith     nrecvs--;
525f32d5d43SBarry Smith     /* unpack a message at a time */
526f32d5d43SBarry Smith     for (j=0; j<rstarts[imdex+1]-rstarts[imdex]; j++) {
527f32d5d43SBarry Smith       for (k=0; k<ncols; k++) {
528f32d5d43SBarry Smith         w[rindices[rstarts[imdex]+j] + nrows*k] = rvalues[ncols*(rstarts[imdex] + j) + k];
5292636ff26SBarry Smith       }
5302636ff26SBarry Smith     }
531f32d5d43SBarry Smith   }
532cb9801acSJed Brown   if (to->n) {ierr = MPI_Waitall(to->n,swaits,to->sstatus);CHKERRQ(ierr);}
533f32d5d43SBarry Smith 
5348c778c55SBarry Smith   ierr = MatDenseRestoreArray(B,&b);CHKERRQ(ierr);
5358c778c55SBarry Smith   ierr = MatDenseRestoreArray(workB,&w);CHKERRQ(ierr);
536f32d5d43SBarry Smith   ierr = MatAssemblyBegin(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
537f32d5d43SBarry Smith   ierr = MatAssemblyEnd(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
538f32d5d43SBarry Smith   PetscFunctionReturn(0);
539f32d5d43SBarry Smith }
540f32d5d43SBarry Smith extern PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat,Mat,Mat);
541f32d5d43SBarry Smith 
5429123193aSHong Zhang #undef __FUNCT__
5439123193aSHong Zhang #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIDense"
5449123193aSHong Zhang PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIDense(Mat A,Mat B,Mat C)
5459123193aSHong Zhang {
5469123193aSHong Zhang   PetscErrorCode ierr;
547f32d5d43SBarry Smith   Mat_MPIAIJ     *aij    = (Mat_MPIAIJ*)A->data;
548f32d5d43SBarry Smith   Mat_MPIDense   *bdense = (Mat_MPIDense*)B->data;
549f32d5d43SBarry Smith   Mat_MPIDense   *cdense = (Mat_MPIDense*)C->data;
550f32d5d43SBarry Smith   Mat            workB;
5519123193aSHong Zhang 
5529123193aSHong Zhang   PetscFunctionBegin;
553f32d5d43SBarry Smith   /* diagonal block of A times all local rows of B*/
554f32d5d43SBarry Smith   ierr = MatMatMultNumeric_SeqAIJ_SeqDense(aij->A,bdense->A,cdense->A);CHKERRQ(ierr);
555f32d5d43SBarry Smith 
556f32d5d43SBarry Smith   /* get off processor parts of B needed to complete the product */
5574324174eSBarry Smith   ierr = MatMPIDenseScatter(A,B,C,&workB);CHKERRQ(ierr);
558f32d5d43SBarry Smith 
559f32d5d43SBarry Smith   /* off-diagonal block of A times nonlocal rows of B */
560f32d5d43SBarry Smith   ierr = MatMatMultNumericAdd_SeqAIJ_SeqDense(aij->B,workB,cdense->A);CHKERRQ(ierr);
5619123193aSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
5629123193aSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
5639123193aSHong Zhang   PetscFunctionReturn(0);
5649123193aSHong Zhang }
565cf1a0672SHong Zhang 
5661634b032SHong Zhang #undef __FUNCT__
567f8efd71cSHong Zhang #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIAIJ"
568f8efd71cSHong Zhang PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C)
5691634b032SHong Zhang {
5701634b032SHong Zhang   PetscErrorCode ierr;
571cf1a0672SHong Zhang   Mat_MPIAIJ     *a   = (Mat_MPIAIJ*)A->data,*c=(Mat_MPIAIJ*)C->data;
572cf1a0672SHong Zhang   Mat_SeqAIJ     *ad  = (Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
573cf1a0672SHong Zhang   Mat_SeqAIJ     *cd  = (Mat_SeqAIJ*)(c->A)->data,*co=(Mat_SeqAIJ*)(c->B)->data;
574cf1a0672SHong Zhang   PetscInt       *adi = ad->i,*adj,*aoi=ao->i,*aoj;
575cf1a0672SHong Zhang   PetscScalar    *ada,*aoa,*cda=cd->a,*coa=co->a;
576cf1a0672SHong Zhang   Mat_SeqAIJ     *p_loc,*p_oth;
577cf1a0672SHong Zhang   PetscInt       *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pj;
57829825010SHong Zhang   PetscScalar    *pa_loc,*pa_oth,*pa,valtmp,*ca;
57929825010SHong Zhang   PetscInt       cm          = C->rmap->n,anz,pnz;
580cf1a0672SHong Zhang   Mat_PtAPMPI    *ptap       = c->ptap;
58129825010SHong Zhang   PetscScalar    *apa_sparse = ptap->apa;
582cf1a0672SHong Zhang   PetscInt       *api,*apj,*apJ,i,j,k,row;
58329825010SHong Zhang   PetscInt       cstart = C->cmap->rstart;
58429825010SHong Zhang   PetscInt       cdnz,conz,k0,k1,nextp;
585a7c7454dSHong Zhang   MPI_Comm       comm;
586a7c7454dSHong Zhang   PetscMPIInt    size;
5871634b032SHong Zhang 
5881634b032SHong Zhang   PetscFunctionBegin;
589a7c7454dSHong Zhang   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
590a7c7454dSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
591a7c7454dSHong Zhang 
592cf1a0672SHong Zhang   /* 1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
593cf1a0672SHong Zhang   /*-----------------------------------------------------*/
594cf1a0672SHong Zhang   /* update numerical values of P_oth and P_loc */
595b7f45c76SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
596cf1a0672SHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
597cf1a0672SHong Zhang 
598cf1a0672SHong Zhang   /* 2) compute numeric C_loc = A_loc*P = Ad*P_loc + Ao*P_oth */
599cf1a0672SHong Zhang   /*----------------------------------------------------------*/
600cf1a0672SHong Zhang   /* get data from symbolic products */
601cf1a0672SHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
602cf1a0672SHong Zhang   pi_loc = p_loc->i; pj_loc = p_loc->j; pa_loc = p_loc->a;
603a7c7454dSHong Zhang   if (size >1) {
604a7c7454dSHong Zhang     p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
605cf1a0672SHong Zhang     pi_oth = p_oth->i; pj_oth = p_oth->j; pa_oth = p_oth->a;
60620e3dc75SHong Zhang   } else {
60720e3dc75SHong Zhang     p_oth = NULL; pi_oth = NULL; pj_oth = NULL; pa_oth = NULL;
608a7c7454dSHong Zhang   }
609cf1a0672SHong Zhang 
61029825010SHong Zhang   api = ptap->api;
61129825010SHong Zhang   apj = ptap->apj;
612cf1a0672SHong Zhang   for (i=0; i<cm; i++) {
61329825010SHong Zhang     apJ = apj + api[i];
61429825010SHong Zhang 
615cf1a0672SHong Zhang     /* diagonal portion of A */
616cf1a0672SHong Zhang     anz = adi[i+1] - adi[i];
617cf1a0672SHong Zhang     adj = ad->j + adi[i];
618cf1a0672SHong Zhang     ada = ad->a + adi[i];
619cf1a0672SHong Zhang     for (j=0; j<anz; j++) {
620cf1a0672SHong Zhang       row = adj[j];
621cf1a0672SHong Zhang       pnz = pi_loc[row+1] - pi_loc[row];
622cf1a0672SHong Zhang       pj  = pj_loc + pi_loc[row];
623cf1a0672SHong Zhang       pa  = pa_loc + pi_loc[row];
62429825010SHong Zhang       /* perform sparse axpy */
625cf1a0672SHong Zhang       valtmp = ada[j];
62629825010SHong Zhang       nextp  = 0;
62729825010SHong Zhang       for (k=0; nextp<pnz; k++) {
62829825010SHong Zhang         if (apJ[k] == pj[nextp]) { /* column of AP == column of P */
62929825010SHong Zhang           apa_sparse[k] += valtmp*pa[nextp++];
63029825010SHong Zhang         }
6311634b032SHong Zhang       }
632cf1a0672SHong Zhang       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
633cf1a0672SHong Zhang     }
6341634b032SHong Zhang 
635cf1a0672SHong Zhang     /* off-diagonal portion of A */
636cf1a0672SHong Zhang     anz = aoi[i+1] - aoi[i];
637cf1a0672SHong Zhang     aoj = ao->j + aoi[i];
638cf1a0672SHong Zhang     aoa = ao->a + aoi[i];
639cf1a0672SHong Zhang     for (j=0; j<anz; j++) {
640cf1a0672SHong Zhang       row = aoj[j];
641cf1a0672SHong Zhang       pnz = pi_oth[row+1] - pi_oth[row];
642cf1a0672SHong Zhang       pj  = pj_oth + pi_oth[row];
643cf1a0672SHong Zhang       pa  = pa_oth + pi_oth[row];
64429825010SHong Zhang       /* perform sparse axpy */
645cf1a0672SHong Zhang       valtmp = aoa[j];
64629825010SHong Zhang       nextp  = 0;
64729825010SHong Zhang       for (k=0; nextp<pnz; k++) {
64829825010SHong Zhang         if (apJ[k] == pj[nextp]) { /* column of AP == column of P */
64929825010SHong Zhang           apa_sparse[k] += valtmp*pa[nextp++];
65029825010SHong Zhang         }
651cf1a0672SHong Zhang       }
652cf1a0672SHong Zhang       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
653cf1a0672SHong Zhang     }
6541634b032SHong Zhang 
655cf1a0672SHong Zhang     /* set values in C */
656cf1a0672SHong Zhang     cdnz = cd->i[i+1] - cd->i[i];
657cf1a0672SHong Zhang     conz = co->i[i+1] - co->i[i];
6581634b032SHong Zhang 
659cf1a0672SHong Zhang     /* 1st off-diagoanl part of C */
660cf1a0672SHong Zhang     ca = coa + co->i[i];
661cf1a0672SHong Zhang     k  = 0;
662cf1a0672SHong Zhang     for (k0=0; k0<conz; k0++) {
663cf1a0672SHong Zhang       if (apJ[k] >= cstart) break;
66429825010SHong Zhang       ca[k0]        = apa_sparse[k];
66529825010SHong Zhang       apa_sparse[k] = 0.0;
666cf1a0672SHong Zhang       k++;
667cf1a0672SHong Zhang     }
6681634b032SHong Zhang 
669cf1a0672SHong Zhang     /* diagonal part of C */
670cf1a0672SHong Zhang     ca = cda + cd->i[i];
671cf1a0672SHong Zhang     for (k1=0; k1<cdnz; k1++) {
67229825010SHong Zhang       ca[k1]        = apa_sparse[k];
67329825010SHong Zhang       apa_sparse[k] = 0.0;
674cf1a0672SHong Zhang       k++;
675cf1a0672SHong Zhang     }
676cf1a0672SHong Zhang 
677cf1a0672SHong Zhang     /* 2nd off-diagoanl part of C */
678cf1a0672SHong Zhang     ca = coa + co->i[i];
679cf1a0672SHong Zhang     for (; k0<conz; k0++) {
68029825010SHong Zhang       ca[k0]        = apa_sparse[k];
68129825010SHong Zhang       apa_sparse[k] = 0.0;
682cf1a0672SHong Zhang       k++;
683cf1a0672SHong Zhang     }
684cf1a0672SHong Zhang   }
685cf1a0672SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
686cf1a0672SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
687cf1a0672SHong Zhang   PetscFunctionReturn(0);
688cf1a0672SHong Zhang }
689cf1a0672SHong Zhang 
6900fc8cf34SHong Zhang /* same as MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(), except using LLCondensed to avoid O(BN) memory requirement */
691cf1a0672SHong Zhang #undef __FUNCT__
692b2405163SHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ"
693b2405163SHong Zhang PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
694cf1a0672SHong Zhang {
695cf1a0672SHong Zhang   PetscErrorCode     ierr;
696ce94432eSBarry Smith   MPI_Comm           comm;
697a7c7454dSHong Zhang   PetscMPIInt        size;
698cf1a0672SHong Zhang   Mat                Cmpi;
699cf1a0672SHong Zhang   Mat_PtAPMPI        *ptap;
7000298fd71SBarry Smith   PetscFreeSpaceList free_space = NULL,current_space=NULL;
701cf1a0672SHong Zhang   Mat_MPIAIJ         *a         = (Mat_MPIAIJ*)A->data,*c;
702cf1a0672SHong Zhang   Mat_SeqAIJ         *ad        = (Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_loc,*p_oth;
703cf1a0672SHong Zhang   PetscInt           *pi_loc,*pj_loc,*pi_oth,*pj_oth,*dnz,*onz;
704cf1a0672SHong Zhang   PetscInt           *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,rstart=A->rmap->rstart;
7050ca7d551SHong Zhang   PetscInt           i,pnz,row,*api,*apj,*Jptr,apnz,nspacedouble=0,j,nzi,*lnk,apnz_max;
7060ca7d551SHong Zhang   PetscInt           am=A->rmap->n,pN=P->cmap->N,pn=P->cmap->n,pm=P->rmap->n;
707cf1a0672SHong Zhang   PetscReal          afill;
708cf1a0672SHong Zhang   PetscScalar        *apa;
709eca6b86aSHong Zhang   PetscTable         ta;
710cf1a0672SHong Zhang 
711cf1a0672SHong Zhang   PetscFunctionBegin;
712ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
713a7c7454dSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
714a7c7454dSHong Zhang 
715cf1a0672SHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
716b00a9115SJed Brown   ierr = PetscNew(&ptap);CHKERRQ(ierr);
717cf1a0672SHong Zhang 
718cf1a0672SHong Zhang   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
719b7f45c76SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
72019f0e57aSHong Zhang 
721cf1a0672SHong Zhang   /* get P_loc by taking all local rows of P */
722cf1a0672SHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
723cf1a0672SHong Zhang 
724cf1a0672SHong Zhang   p_loc  = (Mat_SeqAIJ*)(ptap->P_loc)->data;
725cf1a0672SHong Zhang   pi_loc = p_loc->i; pj_loc = p_loc->j;
726a7c7454dSHong Zhang   if (size > 1) {
727a7c7454dSHong Zhang     p_oth  = (Mat_SeqAIJ*)(ptap->P_oth)->data;
72820e3dc75SHong Zhang     pi_oth = p_oth->i; pj_oth = p_oth->j;
729968382fdSHong Zhang   } else {
73020e3dc75SHong Zhang     p_oth  = NULL;
73120e3dc75SHong Zhang     pi_oth = NULL; pj_oth = NULL;
732a7c7454dSHong Zhang   }
733cf1a0672SHong Zhang 
734cf1a0672SHong Zhang   /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */
735cf1a0672SHong Zhang   /*-------------------------------------------------------------------*/
736854ce69bSBarry Smith   ierr      = PetscMalloc1(am+2,&api);CHKERRQ(ierr);
737cf1a0672SHong Zhang   ptap->api = api;
738cf1a0672SHong Zhang   api[0]    = 0;
739cf1a0672SHong Zhang 
740cf1a0672SHong Zhang   /* create and initialize a linked list */
74145d00d1dSHong Zhang   ierr = PetscTableCreate(pn,pN,&ta);CHKERRQ(ierr);
7420ca7d551SHong Zhang 
7430ca7d551SHong Zhang   /* Calculate apnz_max */
7440ca7d551SHong Zhang   apnz_max = 0;
7450ca7d551SHong Zhang   for (i=0; i<am; i++) {
7460ca7d551SHong Zhang     ierr = PetscTableRemoveAll(ta);CHKERRQ(ierr);
7470ca7d551SHong Zhang     /* diagonal portion of A */
7480ca7d551SHong Zhang     nzi  = adi[i+1] - adi[i];
7490ca7d551SHong Zhang     Jptr = adj+adi[i];  /* cols of A_diag */
7500ca7d551SHong Zhang     MatMergeRows_SeqAIJ(p_loc,nzi,Jptr,ta);
7510ca7d551SHong Zhang     ierr = PetscTableGetCount(ta,&apnz);CHKERRQ(ierr);
7520ca7d551SHong Zhang     if (apnz_max < apnz) apnz_max = apnz;
7530ca7d551SHong Zhang 
7540ca7d551SHong Zhang     /*  off-diagonal portion of A */
7550ca7d551SHong Zhang     nzi = aoi[i+1] - aoi[i];
7560ca7d551SHong Zhang     Jptr = aoj+aoi[i];  /* cols of A_off */
7570ca7d551SHong Zhang     MatMergeRows_SeqAIJ(p_oth,nzi,Jptr,ta);
7580ca7d551SHong Zhang     ierr = PetscTableGetCount(ta,&apnz);CHKERRQ(ierr);
7590ca7d551SHong Zhang     if (apnz_max < apnz) apnz_max = apnz;
7600ca7d551SHong Zhang   }
761eca6b86aSHong Zhang   ierr = PetscTableDestroy(&ta);CHKERRQ(ierr);
762eca6b86aSHong Zhang 
7630ca7d551SHong Zhang   ierr = PetscLLCondensedCreate_Scalable(apnz_max,&lnk);CHKERRQ(ierr);
764cf1a0672SHong Zhang 
765cf1a0672SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */
766f91af8c7SBarry Smith   ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(adi[am],PetscIntSumTruncate(aoi[am],pi_loc[pm]))),&free_space);CHKERRQ(ierr);
767cf1a0672SHong Zhang   current_space = free_space;
768cf1a0672SHong Zhang   ierr = MatPreallocateInitialize(comm,am,pn,dnz,onz);CHKERRQ(ierr);
769cf1a0672SHong Zhang   for (i=0; i<am; i++) {
770cf1a0672SHong Zhang     /* diagonal portion of A */
771cf1a0672SHong Zhang     nzi = adi[i+1] - adi[i];
772cf1a0672SHong Zhang     for (j=0; j<nzi; j++) {
773cf1a0672SHong Zhang       row  = *adj++;
774cf1a0672SHong Zhang       pnz  = pi_loc[row+1] - pi_loc[row];
775cf1a0672SHong Zhang       Jptr = pj_loc + pi_loc[row];
776cf1a0672SHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
777f84c536bSHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr);
778cf1a0672SHong Zhang     }
779cf1a0672SHong Zhang     /* off-diagonal portion of A */
780cf1a0672SHong Zhang     nzi = aoi[i+1] - aoi[i];
781cf1a0672SHong Zhang     for (j=0; j<nzi; j++) {
782cf1a0672SHong Zhang       row  = *aoj++;
783cf1a0672SHong Zhang       pnz  = pi_oth[row+1] - pi_oth[row];
784cf1a0672SHong Zhang       Jptr = pj_oth + pi_oth[row];
785f84c536bSHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr);
786cf1a0672SHong Zhang     }
787cf1a0672SHong Zhang 
788f84c536bSHong Zhang     apnz     = *lnk;
789cf1a0672SHong Zhang     api[i+1] = api[i] + apnz;
790cf1a0672SHong Zhang 
791cf1a0672SHong Zhang     /* if free space is not available, double the total space in the list */
792cf1a0672SHong Zhang     if (current_space->local_remaining<apnz) {
793f91af8c7SBarry Smith       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
794cf1a0672SHong Zhang       nspacedouble++;
795cf1a0672SHong Zhang     }
796cf1a0672SHong Zhang 
797cf1a0672SHong Zhang     /* Copy data into free space, then initialize lnk */
798f84c536bSHong Zhang     ierr = PetscLLCondensedClean_Scalable(apnz,current_space->array,lnk);CHKERRQ(ierr);
799cf1a0672SHong Zhang     ierr = MatPreallocateSet(i+rstart,apnz,current_space->array,dnz,onz);CHKERRQ(ierr);
8002205254eSKarl Rupp 
801cf1a0672SHong Zhang     current_space->array           += apnz;
802cf1a0672SHong Zhang     current_space->local_used      += apnz;
803cf1a0672SHong Zhang     current_space->local_remaining -= apnz;
804cf1a0672SHong Zhang   }
805cf1a0672SHong Zhang 
806cf1a0672SHong Zhang   /* Allocate space for apj, initialize apj, and */
807cf1a0672SHong Zhang   /* destroy list of free space and other temporary array(s) */
808854ce69bSBarry Smith   ierr = PetscMalloc1(api[am]+1,&ptap->apj);CHKERRQ(ierr);
809cf1a0672SHong Zhang   apj  = ptap->apj;
810cf1a0672SHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,ptap->apj);CHKERRQ(ierr);
811f84c536bSHong Zhang   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr);
812cf1a0672SHong Zhang 
813cf1a0672SHong Zhang   /* create and assemble symbolic parallel matrix Cmpi */
814cf1a0672SHong Zhang   /*----------------------------------------------------*/
815cf1a0672SHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
816cf1a0672SHong Zhang   ierr = MatSetSizes(Cmpi,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
81733d57670SJed Brown   ierr = MatSetBlockSizesFromMats(Cmpi,A,P);CHKERRQ(ierr);
818cf1a0672SHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
819cf1a0672SHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
820cf1a0672SHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
821cf1a0672SHong Zhang 
82229825010SHong Zhang   /* malloc apa for assembly Cmpi */
8231795a4d1SJed Brown   ierr = PetscCalloc1(apnz_max,&apa);CHKERRQ(ierr);
8242205254eSKarl Rupp 
82529825010SHong Zhang   ptap->apa = apa;
826cf1a0672SHong Zhang   for (i=0; i<am; i++) {
827cf1a0672SHong Zhang     row  = i + rstart;
828cf1a0672SHong Zhang     apnz = api[i+1] - api[i];
829cf1a0672SHong Zhang     ierr = MatSetValues(Cmpi,1,&row,apnz,apj,apa,INSERT_VALUES);CHKERRQ(ierr);
830cf1a0672SHong Zhang     apj += apnz;
831cf1a0672SHong Zhang   }
832cf1a0672SHong Zhang   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
833cf1a0672SHong Zhang   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
834cf1a0672SHong Zhang 
835cf1a0672SHong Zhang   ptap->destroy             = Cmpi->ops->destroy;
836cf1a0672SHong Zhang   ptap->duplicate           = Cmpi->ops->duplicate;
837f8efd71cSHong Zhang   Cmpi->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ;
838cf1a0672SHong Zhang   Cmpi->ops->destroy        = MatDestroy_MPIAIJ_MatMatMult;
839cf1a0672SHong Zhang   Cmpi->ops->duplicate      = MatDuplicate_MPIAIJ_MatMatMult;
840cf1a0672SHong Zhang 
841cf1a0672SHong Zhang   /* attach the supporting struct to Cmpi for reuse */
842cf1a0672SHong Zhang   c       = (Mat_MPIAIJ*)Cmpi->data;
843cf1a0672SHong Zhang   c->ptap = ptap;
844cf1a0672SHong Zhang 
845cf1a0672SHong Zhang   *C = Cmpi;
846cf1a0672SHong Zhang 
847cf1a0672SHong Zhang   /* set MatInfo */
848118727c9SMark F. Adams   afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1) + 1.e-5;
849cf1a0672SHong Zhang   if (afill < 1.0) afill = 1.0;
850cf1a0672SHong Zhang   Cmpi->info.mallocs           = nspacedouble;
851cf1a0672SHong Zhang   Cmpi->info.fill_ratio_given  = fill;
852cf1a0672SHong Zhang   Cmpi->info.fill_ratio_needed = afill;
853cf1a0672SHong Zhang 
854cf1a0672SHong Zhang #if defined(PETSC_USE_INFO)
855cf1a0672SHong Zhang   if (api[am]) {
85657622a8eSBarry Smith     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
85757622a8eSBarry Smith     ierr = PetscInfo1(Cmpi,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr);
858cf1a0672SHong Zhang   } else {
859cf1a0672SHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
860cf1a0672SHong Zhang   }
861cf1a0672SHong Zhang #endif
8621634b032SHong Zhang   PetscFunctionReturn(0);
8631634b032SHong Zhang }
864f96d37fcSHong Zhang 
865f96d37fcSHong Zhang /*-------------------------------------------------------------------------*/
866f96d37fcSHong Zhang #undef __FUNCT__
867f96d37fcSHong Zhang #define __FUNCT__ "MatTransposeMatMult_MPIAIJ_MPIAIJ"
868c5af039cSHong Zhang PetscErrorCode MatTransposeMatMult_MPIAIJ_MPIAIJ(Mat P,Mat A,MatReuse scall,PetscReal fill,Mat *C)
869f96d37fcSHong Zhang {
870f96d37fcSHong Zhang   PetscErrorCode ierr;
871c216dbf3SHong Zhang   const char     *algTypes[3] = {"scalable","nonscalable","matmatmult"};
872c216dbf3SHong Zhang   PetscInt       alg=0; /* set default algorithm */
873f96d37fcSHong Zhang 
874f96d37fcSHong Zhang   PetscFunctionBegin;
875c216dbf3SHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
876c216dbf3SHong Zhang     ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr);
877c216dbf3SHong Zhang     ierr = PetscOptionsEList("-mattransposematmult_via","Algorithmic approach","MatTransposeMatMult",algTypes,3,algTypes[0],&alg,NULL);CHKERRQ(ierr);
878c216dbf3SHong Zhang     ierr = PetscOptionsEnd();CHKERRQ(ierr);
879c216dbf3SHong Zhang 
880c216dbf3SHong Zhang     ierr = PetscLogEventBegin(MAT_TransposeMatMultSymbolic,P,A,0,0);CHKERRQ(ierr);
881c216dbf3SHong Zhang     switch (alg) {
882c216dbf3SHong Zhang     case 1:
8832bbb1c24SHong Zhang       ierr = MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(P,A,fill,C);CHKERRQ(ierr);
884c216dbf3SHong Zhang       break;
885c216dbf3SHong Zhang     case 2:
886275476c6SMatthew G. Knepley     {
887ab1b013aSHong Zhang       Mat         Pt;
888ab1b013aSHong Zhang       Mat_PtAPMPI *ptap;
889ab1b013aSHong Zhang       Mat_MPIAIJ  *c;
890ab1b013aSHong Zhang       ierr = MatTranspose(P,MAT_INITIAL_MATRIX,&Pt);CHKERRQ(ierr);
891ab1b013aSHong Zhang       ierr = MatMatMult(Pt,A,MAT_INITIAL_MATRIX,fill,C);CHKERRQ(ierr);
892ab1b013aSHong Zhang       c        = (Mat_MPIAIJ*)(*C)->data;
893ab1b013aSHong Zhang       ptap     = c->ptap;
894ab1b013aSHong Zhang       ptap->Pt = Pt;
895c216dbf3SHong Zhang       (*C)->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult;
896c216dbf3SHong Zhang       PetscFunctionReturn(0);
897275476c6SMatthew G. Knepley     }
898c216dbf3SHong Zhang       break;
899c216dbf3SHong Zhang     default:
9006da69ca6SHong Zhang       ierr = MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(P,A,fill,C);CHKERRQ(ierr);
901c216dbf3SHong Zhang       break;
902c216dbf3SHong Zhang     }
903c216dbf3SHong Zhang     ierr = PetscLogEventEnd(MAT_TransposeMatMultSymbolic,P,A,0,0);CHKERRQ(ierr);
904c216dbf3SHong Zhang   }
905c216dbf3SHong Zhang   ierr = PetscLogEventBegin(MAT_TransposeMatMultNumeric,P,A,0,0);CHKERRQ(ierr);
906c216dbf3SHong Zhang   ierr = (*(*C)->ops->mattransposemultnumeric)(P,A,*C);CHKERRQ(ierr);
907c216dbf3SHong Zhang   ierr = PetscLogEventEnd(MAT_TransposeMatMultNumeric,P,A,0,0);CHKERRQ(ierr);
908ab1b013aSHong Zhang   PetscFunctionReturn(0);
909ab1b013aSHong Zhang }
910ab1b013aSHong Zhang 
911c216dbf3SHong Zhang /* This routine only works when scall=MAT_REUSE_MATRIX! */
912c216dbf3SHong Zhang #undef __FUNCT__
913c216dbf3SHong Zhang #define __FUNCT__ "MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult"
914c216dbf3SHong Zhang PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult(Mat P,Mat A,Mat C)
915c216dbf3SHong Zhang {
916c216dbf3SHong Zhang   PetscErrorCode ierr;
9172bbb1c24SHong Zhang   Mat_MPIAIJ     *c=(Mat_MPIAIJ*)C->data;
9182bbb1c24SHong Zhang   Mat_PtAPMPI    *ptap= c->ptap;
9192bbb1c24SHong Zhang   Mat            Pt=ptap->Pt;
920c216dbf3SHong Zhang 
921c216dbf3SHong Zhang   PetscFunctionBegin;
922c216dbf3SHong Zhang   ierr = MatTranspose(P,MAT_REUSE_MATRIX,&Pt);CHKERRQ(ierr);
923c216dbf3SHong Zhang   ierr = MatMatMultNumeric(Pt,A,C);CHKERRQ(ierr);
924f96d37fcSHong Zhang   PetscFunctionReturn(0);
925f96d37fcSHong Zhang }
926f96d37fcSHong Zhang 
9276da69ca6SHong Zhang /* Non-scalable version, use dense axpy */
928f96d37fcSHong Zhang #undef __FUNCT__
9296da69ca6SHong Zhang #define __FUNCT__ "MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable"
9306da69ca6SHong Zhang PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable(Mat P,Mat A,Mat C)
931f96d37fcSHong Zhang {
932c5af039cSHong Zhang   PetscErrorCode      ierr;
933c5af039cSHong Zhang   Mat_Merge_SeqsToMPI *merge;
934497f5370SHong Zhang   Mat_MPIAIJ          *p =(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
935c5af039cSHong Zhang   Mat_SeqAIJ          *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data;
936c5af039cSHong Zhang   Mat_PtAPMPI         *ptap;
937d6ab1dc0SHong Zhang   PetscInt            *adj,*aJ;
938497f5370SHong Zhang   PetscInt            i,j,k,anz,pnz,row,*cj;
939e745005bSHong Zhang   MatScalar           *ada,*aval,*ca,valtmp;
940c5af039cSHong Zhang   PetscInt            am  =A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n;
941ce94432eSBarry Smith   MPI_Comm            comm;
942c5af039cSHong Zhang   PetscMPIInt         size,rank,taga,*len_s;
943c5af039cSHong Zhang   PetscInt            *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci;
944c5af039cSHong Zhang   PetscInt            **buf_ri,**buf_rj;
945c5af039cSHong Zhang   PetscInt            cnz=0,*bj_i,*bi,*bj,bnz,nextcj;  /* bi,bj,ba: local array of C(mpi mat) */
946c5af039cSHong Zhang   MPI_Request         *s_waits,*r_waits;
947c5af039cSHong Zhang   MPI_Status          *status;
948c5af039cSHong Zhang   MatScalar           **abuf_r,*ba_i,*pA,*coa,*ba;
949d6ab1dc0SHong Zhang   PetscInt            *ai,*aj,*coi,*coj;
950a2ea699eSBarry Smith   PetscInt            *poJ,*pdJ;
951c5af039cSHong Zhang   Mat                 A_loc;
952c5af039cSHong Zhang   Mat_SeqAIJ          *a_loc;
953f96d37fcSHong Zhang 
954f96d37fcSHong Zhang   PetscFunctionBegin;
955ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr);
956c5af039cSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
957c5af039cSHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
958c5af039cSHong Zhang 
959c5af039cSHong Zhang   ptap  = c->ptap;
960c5af039cSHong Zhang   merge = ptap->merge;
961c5af039cSHong Zhang 
962c5af039cSHong Zhang   /* 2) compute numeric C_seq = P_loc^T*A_loc*P - dominating part */
963c5af039cSHong Zhang   /*--------------------------------------------------------------*/
964c5af039cSHong Zhang   /* get data from symbolic products */
965c5af039cSHong Zhang   coi  = merge->coi; coj = merge->coj;
966854ce69bSBarry Smith   ierr = PetscCalloc1(coi[pon]+1,&coa);CHKERRQ(ierr);
967c5af039cSHong Zhang 
968c5af039cSHong Zhang   bi     = merge->bi; bj = merge->bj;
969c5af039cSHong Zhang   owners = merge->rowmap->range;
970854ce69bSBarry Smith   ierr   = PetscCalloc1(bi[cm]+1,&ba);CHKERRQ(ierr);
971c5af039cSHong Zhang 
972c5af039cSHong Zhang   /* get A_loc by taking all local rows of A */
973c5af039cSHong Zhang   A_loc = ptap->A_loc;
974c5af039cSHong Zhang   ierr  = MatMPIAIJGetLocalMat(A,MAT_REUSE_MATRIX,&A_loc);CHKERRQ(ierr);
975c5af039cSHong Zhang   a_loc = (Mat_SeqAIJ*)(A_loc)->data;
976d6ab1dc0SHong Zhang   ai    = a_loc->i;
977d6ab1dc0SHong Zhang   aj    = a_loc->j;
978c5af039cSHong Zhang 
979854ce69bSBarry Smith   ierr = PetscCalloc1(A->cmap->N,&aval);CHKERRQ(ierr); /* non-scalable!!! */
980c5af039cSHong Zhang 
981c5af039cSHong Zhang   for (i=0; i<am; i++) {
982e745005bSHong Zhang     /* 2-a) put A[i,:] to dense array aval */
983d6ab1dc0SHong Zhang     anz = ai[i+1] - ai[i];
984d6ab1dc0SHong Zhang     adj = aj + ai[i];
985d6ab1dc0SHong Zhang     ada = a_loc->a + ai[i];
986c5af039cSHong Zhang     for (j=0; j<anz; j++) {
987e745005bSHong Zhang       aval[adj[j]] = ada[j];
988c5af039cSHong Zhang     }
989c5af039cSHong Zhang 
990c5af039cSHong Zhang     /* 2-b) Compute Cseq = P_loc[i,:]^T*A[i,:] using outer product */
991c5af039cSHong Zhang     /*--------------------------------------------------------------*/
992c5af039cSHong Zhang     /* put the value into Co=(p->B)^T*A (off-diagonal part, send to others) */
993c5af039cSHong Zhang     pnz = po->i[i+1] - po->i[i];
994c5af039cSHong Zhang     poJ = po->j + po->i[i];
995c5af039cSHong Zhang     pA  = po->a + po->i[i];
996c5af039cSHong Zhang     for (j=0; j<pnz; j++) {
997c5af039cSHong Zhang       row = poJ[j];
998c5af039cSHong Zhang       cnz = coi[row+1] - coi[row];
999c5af039cSHong Zhang       cj  = coj + coi[row];
1000c5af039cSHong Zhang       ca  = coa + coi[row];
1001c5af039cSHong Zhang       /* perform dense axpy */
1002c5af039cSHong Zhang       valtmp = pA[j];
1003c5af039cSHong Zhang       for (k=0; k<cnz; k++) {
1004e745005bSHong Zhang         ca[k] += valtmp*aval[cj[k]];
1005c5af039cSHong Zhang       }
1006c5af039cSHong Zhang       ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1007c5af039cSHong Zhang     }
1008c5af039cSHong Zhang 
1009c5af039cSHong Zhang     /* put the value into Cd (diagonal part) */
1010c5af039cSHong Zhang     pnz = pd->i[i+1] - pd->i[i];
1011c5af039cSHong Zhang     pdJ = pd->j + pd->i[i];
1012c5af039cSHong Zhang     pA  = pd->a + pd->i[i];
1013c5af039cSHong Zhang     for (j=0; j<pnz; j++) {
1014c5af039cSHong Zhang       row = pdJ[j];
1015c5af039cSHong Zhang       cnz = bi[row+1] - bi[row];
1016c5af039cSHong Zhang       cj  = bj + bi[row];
1017c5af039cSHong Zhang       ca  = ba + bi[row];
1018c5af039cSHong Zhang       /* perform dense axpy */
1019c5af039cSHong Zhang       valtmp = pA[j];
1020c5af039cSHong Zhang       for (k=0; k<cnz; k++) {
1021e745005bSHong Zhang         ca[k] += valtmp*aval[cj[k]];
1022c5af039cSHong Zhang       }
1023c5af039cSHong Zhang       ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1024c5af039cSHong Zhang     }
1025c5af039cSHong Zhang 
1026d6ab1dc0SHong Zhang     /* zero the current row of Pt*A */
1027d6ab1dc0SHong Zhang     aJ = aj + ai[i];
1028e745005bSHong Zhang     for (k=0; k<anz; k++) aval[aJ[k]] = 0.0;
1029c5af039cSHong Zhang   }
1030c5af039cSHong Zhang 
1031c5af039cSHong Zhang   /* 3) send and recv matrix values coa */
1032c5af039cSHong Zhang   /*------------------------------------*/
1033c5af039cSHong Zhang   buf_ri = merge->buf_ri;
1034c5af039cSHong Zhang   buf_rj = merge->buf_rj;
1035c5af039cSHong Zhang   len_s  = merge->len_s;
1036c5af039cSHong Zhang   ierr   = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr);
1037c5af039cSHong Zhang   ierr   = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr);
1038c5af039cSHong Zhang 
1039dcca6d9dSJed Brown   ierr = PetscMalloc2(merge->nsend+1,&s_waits,size,&status);CHKERRQ(ierr);
1040c5af039cSHong Zhang   for (proc=0,k=0; proc<size; proc++) {
1041c5af039cSHong Zhang     if (!len_s[proc]) continue;
1042c5af039cSHong Zhang     i    = merge->owners_co[proc];
1043c5af039cSHong Zhang     ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr);
1044c5af039cSHong Zhang     k++;
1045c5af039cSHong Zhang   }
1046c5af039cSHong Zhang   if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);}
1047c5af039cSHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);}
1048c5af039cSHong Zhang 
1049c5af039cSHong Zhang   ierr = PetscFree2(s_waits,status);CHKERRQ(ierr);
1050c5af039cSHong Zhang   ierr = PetscFree(r_waits);CHKERRQ(ierr);
1051c5af039cSHong Zhang   ierr = PetscFree(coa);CHKERRQ(ierr);
1052c5af039cSHong Zhang 
1053c5af039cSHong Zhang   /* 4) insert local Cseq and received values into Cmpi */
1054c5af039cSHong Zhang   /*----------------------------------------------------*/
1055dcca6d9dSJed Brown   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
1056c5af039cSHong Zhang   for (k=0; k<merge->nrecv; k++) {
1057c36aecf5SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1058c5af039cSHong Zhang     nrows       = *(buf_ri_k[k]);
1059c5af039cSHong Zhang     nextrow[k]  = buf_ri_k[k]+1;  /* next row number of k-th recved i-structure */
1060c5af039cSHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
1061c5af039cSHong Zhang   }
1062c5af039cSHong Zhang 
1063c5af039cSHong Zhang   for (i=0; i<cm; i++) {
1064c5af039cSHong Zhang     row  = owners[rank] + i; /* global row index of C_seq */
1065c5af039cSHong Zhang     bj_i = bj + bi[i];  /* col indices of the i-th row of C */
1066c5af039cSHong Zhang     ba_i = ba + bi[i];
1067c5af039cSHong Zhang     bnz  = bi[i+1] - bi[i];
1068c5af039cSHong Zhang     /* add received vals into ba */
1069c5af039cSHong Zhang     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1070c5af039cSHong Zhang       /* i-th row */
1071c5af039cSHong Zhang       if (i == *nextrow[k]) {
1072c5af039cSHong Zhang         cnz    = *(nextci[k]+1) - *nextci[k];
1073c5af039cSHong Zhang         cj     = buf_rj[k] + *(nextci[k]);
1074c5af039cSHong Zhang         ca     = abuf_r[k] + *(nextci[k]);
1075c5af039cSHong Zhang         nextcj = 0;
1076c5af039cSHong Zhang         for (j=0; nextcj<cnz; j++) {
1077c5af039cSHong Zhang           if (bj_i[j] == cj[nextcj]) { /* bcol == ccol */
1078c5af039cSHong Zhang             ba_i[j] += ca[nextcj++];
1079c5af039cSHong Zhang           }
1080c5af039cSHong Zhang         }
1081c5af039cSHong Zhang         nextrow[k]++; nextci[k]++;
1082c5af039cSHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1083c5af039cSHong Zhang       }
1084c5af039cSHong Zhang     }
1085c5af039cSHong Zhang     ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr);
1086c5af039cSHong Zhang   }
1087c5af039cSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1088c5af039cSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1089c5af039cSHong Zhang 
1090c5af039cSHong Zhang   ierr = PetscFree(ba);CHKERRQ(ierr);
1091c5af039cSHong Zhang   ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr);
1092c5af039cSHong Zhang   ierr = PetscFree(abuf_r);CHKERRQ(ierr);
1093c5af039cSHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1094e745005bSHong Zhang   ierr = PetscFree(aval);CHKERRQ(ierr);
1095f96d37fcSHong Zhang   PetscFunctionReturn(0);
1096f96d37fcSHong Zhang }
1097f96d37fcSHong Zhang 
1098c9ba551fSBarry Smith PetscErrorCode MatDuplicate_MPIAIJ_MatPtAP(Mat, MatDuplicateOption,Mat*);
1099c5af039cSHong Zhang /* This routine is modified from MatPtAPSymbolic_MPIAIJ_MPIAIJ() */
1100f96d37fcSHong Zhang #undef __FUNCT__
11012bbb1c24SHong Zhang #define __FUNCT__ "MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable"
11022bbb1c24SHong Zhang PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat P,Mat A,PetscReal fill,Mat *C)
1103f96d37fcSHong Zhang {
1104f96d37fcSHong Zhang   PetscErrorCode      ierr;
11054e938277SHong Zhang   Mat                 Cmpi,A_loc,POt,PDt;
1106f96d37fcSHong Zhang   Mat_PtAPMPI         *ptap;
11070298fd71SBarry Smith   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
1108f96d37fcSHong Zhang   Mat_MPIAIJ          *p        =(Mat_MPIAIJ*)P->data,*c;
1109f96d37fcSHong Zhang   PetscInt            *pdti,*pdtj,*poti,*potj,*ptJ;
1110f96d37fcSHong Zhang   PetscInt            nnz;
11114e938277SHong Zhang   PetscInt            *lnk,*owners_co,*coi,*coj,i,k,pnz,row;
1112497f5370SHong Zhang   PetscInt            am=A->rmap->n,pn=P->cmap->n;
1113f96d37fcSHong Zhang   PetscBT             lnkbt;
1114ce94432eSBarry Smith   MPI_Comm            comm;
1115f96d37fcSHong Zhang   PetscMPIInt         size,rank,tagi,tagj,*len_si,*len_s,*len_ri;
1116f96d37fcSHong Zhang   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
1117f96d37fcSHong Zhang   PetscInt            len,proc,*dnz,*onz,*owners;
1118f96d37fcSHong Zhang   PetscInt            nzi,*bi,*bj;
1119f96d37fcSHong Zhang   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
1120f96d37fcSHong Zhang   MPI_Request         *swaits,*rwaits;
1121f96d37fcSHong Zhang   MPI_Status          *sstatus,rstatus;
1122f96d37fcSHong Zhang   Mat_Merge_SeqsToMPI *merge;
1123d6ab1dc0SHong Zhang   PetscInt            *ai,*aj,*Jptr,anz,*prmap=p->garray,pon,nspacedouble=0,j;
1124f96d37fcSHong Zhang   PetscReal           afill  =1.0,afill_tmp;
1125438d860cSHong Zhang   PetscInt            rstart = P->cmap->rstart,rmax,aN=A->cmap->N;
1126f96d37fcSHong Zhang   PetscScalar         *vals;
11274e938277SHong Zhang   Mat_SeqAIJ          *a_loc, *pdt,*pot;
1128f96d37fcSHong Zhang 
1129f96d37fcSHong Zhang   PetscFunctionBegin;
1130ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
1131c5af039cSHong Zhang   /* check if matrix local sizes are compatible */
11326c4ed002SBarry Smith   if (A->rmap->rstart != P->rmap->rstart || A->rmap->rend != P->rmap->rend) SETERRQ4(comm,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, A (%D, %D) != P (%D,%D)",A->rmap->rstart,A->rmap->rend,P->rmap->rstart,P->rmap->rend);
1133c5af039cSHong Zhang 
1134f96d37fcSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1135f96d37fcSHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1136f96d37fcSHong Zhang 
1137f96d37fcSHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
1138b00a9115SJed Brown   ierr = PetscNew(&ptap);CHKERRQ(ierr);
1139f96d37fcSHong Zhang 
1140f96d37fcSHong Zhang   /* get A_loc by taking all local rows of A */
1141f96d37fcSHong Zhang   ierr = MatMPIAIJGetLocalMat(A,MAT_INITIAL_MATRIX,&A_loc);CHKERRQ(ierr);
11422205254eSKarl Rupp 
1143c5af039cSHong Zhang   ptap->A_loc = A_loc;
11442205254eSKarl Rupp 
11451c7d5954SHong Zhang   a_loc = (Mat_SeqAIJ*)(A_loc)->data;
1146d6ab1dc0SHong Zhang   ai    = a_loc->i;
1147d6ab1dc0SHong Zhang   aj    = a_loc->j;
1148f96d37fcSHong Zhang 
1149f96d37fcSHong Zhang   /* determine symbolic Co=(p->B)^T*A - send to others */
1150f96d37fcSHong Zhang   /*----------------------------------------------------*/
11514e938277SHong Zhang   ierr = MatTransposeSymbolic_SeqAIJ(p->A,&PDt);CHKERRQ(ierr);
11524e938277SHong Zhang   pdt  = (Mat_SeqAIJ*)PDt->data;
11534e938277SHong Zhang   pdti = pdt->i; pdtj = pdt->j;
11544e938277SHong Zhang 
11554e938277SHong Zhang   ierr = MatTransposeSymbolic_SeqAIJ(p->B,&POt);CHKERRQ(ierr);
11564e938277SHong Zhang   pot  = (Mat_SeqAIJ*)POt->data;
11574e938277SHong Zhang   poti = pot->i; potj = pot->j;
1158f96d37fcSHong Zhang 
1159f96d37fcSHong Zhang   /* then, compute symbolic Co = (p->B)^T*A */
11602205254eSKarl Rupp   pon    = (p->B)->cmap->n; /* total num of rows to be sent to other processors >= (num of nonzero rows of C_seq) - pn */
1161854ce69bSBarry Smith   ierr   = PetscMalloc1(pon+1,&coi);CHKERRQ(ierr);
1162f96d37fcSHong Zhang   coi[0] = 0;
1163f96d37fcSHong Zhang 
1164f96d37fcSHong Zhang   /* set initial free space to be fill*(nnz(p->B) + nnz(A)) */
1165f91af8c7SBarry Smith   nnz           = PetscRealIntMultTruncate(fill,PetscIntSumTruncate(poti[pon],ai[am]));
1166a2ea699eSBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
1167f96d37fcSHong Zhang   current_space = free_space;
116819f0e57aSHong Zhang 
1169f96d37fcSHong Zhang   /* create and initialize a linked list */
1170438d860cSHong Zhang   ierr = PetscLLCondensedCreate(aN,aN,&lnk,&lnkbt);CHKERRQ(ierr);
1171f96d37fcSHong Zhang 
1172f96d37fcSHong Zhang   for (i=0; i<pon; i++) {
1173f96d37fcSHong Zhang     pnz = poti[i+1] - poti[i];
1174f96d37fcSHong Zhang     ptJ = potj + poti[i];
1175f96d37fcSHong Zhang     for (j=0; j<pnz; j++) {
1176f96d37fcSHong Zhang       row  = ptJ[j]; /* row of A_loc == col of Pot */
1177d6ab1dc0SHong Zhang       anz  = ai[row+1] - ai[row];
1178d6ab1dc0SHong Zhang       Jptr = aj + ai[row];
1179f96d37fcSHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
1180d6ab1dc0SHong Zhang       ierr = PetscLLCondensedAddSorted(anz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1181f96d37fcSHong Zhang     }
11824e938277SHong Zhang     nnz = lnk[0];
1183f96d37fcSHong Zhang 
1184f96d37fcSHong Zhang     /* If free space is not available, double the total space in the list */
1185f96d37fcSHong Zhang     if (current_space->local_remaining<nnz) {
1186f91af8c7SBarry Smith       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(nnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
1187f96d37fcSHong Zhang       nspacedouble++;
1188f96d37fcSHong Zhang     }
1189f96d37fcSHong Zhang 
1190f96d37fcSHong Zhang     /* Copy data into free space, and zero out denserows */
11914e938277SHong Zhang     ierr = PetscLLCondensedClean(aN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
11922205254eSKarl Rupp 
1193f96d37fcSHong Zhang     current_space->array           += nnz;
1194f96d37fcSHong Zhang     current_space->local_used      += nnz;
1195f96d37fcSHong Zhang     current_space->local_remaining -= nnz;
11962205254eSKarl Rupp 
1197f96d37fcSHong Zhang     coi[i+1] = coi[i] + nnz;
1198f96d37fcSHong Zhang   }
1199f96d37fcSHong Zhang 
1200854ce69bSBarry Smith   ierr = PetscMalloc1(coi[pon]+1,&coj);CHKERRQ(ierr);
1201f96d37fcSHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr);
12022205254eSKarl Rupp 
1203118727c9SMark F. Adams   afill_tmp = (PetscReal)coi[pon]/(poti[pon] + ai[am]+1);
1204f96d37fcSHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
1205f96d37fcSHong Zhang 
1206f96d37fcSHong Zhang   /* send j-array (coj) of Co to other processors */
1207f96d37fcSHong Zhang   /*----------------------------------------------*/
1208f96d37fcSHong Zhang   /* determine row ownership */
1209b00a9115SJed Brown   ierr = PetscNew(&merge);CHKERRQ(ierr);
1210f96d37fcSHong Zhang   ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr);
12112205254eSKarl Rupp 
1212f96d37fcSHong Zhang   merge->rowmap->n  = pn;
1213f96d37fcSHong Zhang   merge->rowmap->bs = 1;
12142205254eSKarl Rupp 
1215f96d37fcSHong Zhang   ierr   = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr);
1216f96d37fcSHong Zhang   owners = merge->rowmap->range;
1217f96d37fcSHong Zhang 
1218f96d37fcSHong Zhang   /* determine the number of messages to send, their lengths */
12191795a4d1SJed Brown   ierr = PetscCalloc1(size,&len_si);CHKERRQ(ierr);
1220785e854fSJed Brown   ierr = PetscMalloc1(size,&merge->len_s);CHKERRQ(ierr);
12212205254eSKarl Rupp 
1222f96d37fcSHong Zhang   len_s        = merge->len_s;
1223f96d37fcSHong Zhang   merge->nsend = 0;
1224f96d37fcSHong Zhang 
1225854ce69bSBarry Smith   ierr = PetscMalloc1(size+2,&owners_co);CHKERRQ(ierr);
1226f96d37fcSHong Zhang   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
1227f96d37fcSHong Zhang 
1228f96d37fcSHong Zhang   proc = 0;
1229f96d37fcSHong Zhang   for (i=0; i<pon; i++) {
1230f96d37fcSHong Zhang     while (prmap[i] >= owners[proc+1]) proc++;
1231f96d37fcSHong Zhang     len_si[proc]++;  /* num of rows in Co to be sent to [proc] */
1232f96d37fcSHong Zhang     len_s[proc] += coi[i+1] - coi[i];
1233f96d37fcSHong Zhang   }
1234f96d37fcSHong Zhang 
1235f96d37fcSHong Zhang   len          = 0; /* max length of buf_si[] */
1236f96d37fcSHong Zhang   owners_co[0] = 0;
1237f96d37fcSHong Zhang   for (proc=0; proc<size; proc++) {
1238f96d37fcSHong Zhang     owners_co[proc+1] = owners_co[proc] + len_si[proc];
1239f96d37fcSHong Zhang     if (len_si[proc]) {
1240f96d37fcSHong Zhang       merge->nsend++;
1241f96d37fcSHong Zhang       len_si[proc] = 2*(len_si[proc] + 1);
1242f96d37fcSHong Zhang       len         += len_si[proc];
1243f96d37fcSHong Zhang     }
1244f96d37fcSHong Zhang   }
1245f96d37fcSHong Zhang 
1246f96d37fcSHong Zhang   /* determine the number and length of messages to receive for coi and coj  */
12470298fd71SBarry Smith   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&merge->nrecv);CHKERRQ(ierr);
1248f96d37fcSHong Zhang   ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr);
1249f96d37fcSHong Zhang 
1250f96d37fcSHong Zhang   /* post the Irecv and Isend of coj */
1251f96d37fcSHong Zhang   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
1252f96d37fcSHong Zhang   ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
1253854ce69bSBarry Smith   ierr = PetscMalloc1(merge->nsend+1,&swaits);CHKERRQ(ierr);
1254f96d37fcSHong Zhang   for (proc=0, k=0; proc<size; proc++) {
1255f96d37fcSHong Zhang     if (!len_s[proc]) continue;
1256f96d37fcSHong Zhang     i    = owners_co[proc];
1257f96d37fcSHong Zhang     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
1258f96d37fcSHong Zhang     k++;
1259f96d37fcSHong Zhang   }
1260f96d37fcSHong Zhang 
1261f96d37fcSHong Zhang   /* receives and sends of coj are complete */
1262785e854fSJed Brown   ierr = PetscMalloc1(size,&sstatus);CHKERRQ(ierr);
1263f96d37fcSHong Zhang   for (i=0; i<merge->nrecv; i++) {
1264c280ed6aSJed Brown     PetscMPIInt icompleted;
1265c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1266f96d37fcSHong Zhang   }
1267f96d37fcSHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1268f96d37fcSHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1269f96d37fcSHong Zhang 
1270f96d37fcSHong Zhang   /* send and recv coi */
1271f96d37fcSHong Zhang   /*-------------------*/
1272f96d37fcSHong Zhang   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
1273f96d37fcSHong Zhang   ierr   = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
1274854ce69bSBarry Smith   ierr   = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr);
1275f96d37fcSHong Zhang   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
1276f96d37fcSHong Zhang   for (proc=0,k=0; proc<size; proc++) {
1277f96d37fcSHong Zhang     if (!len_s[proc]) continue;
1278f96d37fcSHong Zhang     /* form outgoing message for i-structure:
1279f96d37fcSHong Zhang          buf_si[0]:                 nrows to be sent
1280f96d37fcSHong Zhang                [1:nrows]:           row index (global)
1281f96d37fcSHong Zhang                [nrows+1:2*nrows+1]: i-structure index
1282f96d37fcSHong Zhang     */
1283f96d37fcSHong Zhang     /*-------------------------------------------*/
1284f96d37fcSHong Zhang     nrows       = len_si[proc]/2 - 1;
1285f96d37fcSHong Zhang     buf_si_i    = buf_si + nrows+1;
1286f96d37fcSHong Zhang     buf_si[0]   = nrows;
1287f96d37fcSHong Zhang     buf_si_i[0] = 0;
1288f96d37fcSHong Zhang     nrows       = 0;
1289f96d37fcSHong Zhang     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
1290f96d37fcSHong Zhang       nzi               = coi[i+1] - coi[i];
1291f96d37fcSHong Zhang       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */
1292f96d37fcSHong Zhang       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
1293f96d37fcSHong Zhang       nrows++;
1294f96d37fcSHong Zhang     }
1295f96d37fcSHong Zhang     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
1296f96d37fcSHong Zhang     k++;
1297f96d37fcSHong Zhang     buf_si += len_si[proc];
1298f96d37fcSHong Zhang   }
1299f96d37fcSHong Zhang   i = merge->nrecv;
1300f96d37fcSHong Zhang   while (i--) {
1301c280ed6aSJed Brown     PetscMPIInt icompleted;
1302c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1303f96d37fcSHong Zhang   }
1304f96d37fcSHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1305f96d37fcSHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1306f96d37fcSHong Zhang   ierr = PetscFree(len_si);CHKERRQ(ierr);
1307f96d37fcSHong Zhang   ierr = PetscFree(len_ri);CHKERRQ(ierr);
1308f96d37fcSHong Zhang   ierr = PetscFree(swaits);CHKERRQ(ierr);
1309f96d37fcSHong Zhang   ierr = PetscFree(sstatus);CHKERRQ(ierr);
1310f96d37fcSHong Zhang   ierr = PetscFree(buf_s);CHKERRQ(ierr);
1311f96d37fcSHong Zhang 
1312f96d37fcSHong Zhang   /* compute the local portion of C (mpi mat) */
1313f96d37fcSHong Zhang   /*------------------------------------------*/
1314f96d37fcSHong Zhang   /* allocate bi array and free space for accumulating nonzero column info */
1315854ce69bSBarry Smith   ierr  = PetscMalloc1(pn+1,&bi);CHKERRQ(ierr);
1316f96d37fcSHong Zhang   bi[0] = 0;
1317f96d37fcSHong Zhang 
1318c36aecf5SHong Zhang   /* set initial free space to be fill*(nnz(P) + nnz(A)) */
1319f91af8c7SBarry Smith   nnz           = PetscRealIntMultTruncate(fill,PetscIntSumTruncate(pdti[pn],PetscIntSumTruncate(poti[pon],ai[am])));
1320a2ea699eSBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
1321f96d37fcSHong Zhang   current_space = free_space;
1322f96d37fcSHong Zhang 
1323dcca6d9dSJed Brown   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
1324f96d37fcSHong Zhang   for (k=0; k<merge->nrecv; k++) {
1325f96d37fcSHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1326f96d37fcSHong Zhang     nrows       = *buf_ri_k[k];
1327f96d37fcSHong Zhang     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
1328f96d37fcSHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
1329f96d37fcSHong Zhang   }
1330f96d37fcSHong Zhang 
13311c7d5954SHong Zhang   ierr = MatPreallocateInitialize(comm,pn,A->cmap->n,dnz,onz);CHKERRQ(ierr);
1332f96d37fcSHong Zhang   rmax = 0;
1333f96d37fcSHong Zhang   for (i=0; i<pn; i++) {
1334f96d37fcSHong Zhang     /* add pdt[i,:]*AP into lnk */
1335f96d37fcSHong Zhang     pnz = pdti[i+1] - pdti[i];
1336f96d37fcSHong Zhang     ptJ = pdtj + pdti[i];
1337f96d37fcSHong Zhang     for (j=0; j<pnz; j++) {
1338f96d37fcSHong Zhang       row  = ptJ[j];  /* row of AP == col of Pt */
1339d6ab1dc0SHong Zhang       anz  = ai[row+1] - ai[row];
1340d6ab1dc0SHong Zhang       Jptr = aj + ai[row];
1341f96d37fcSHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
1342d6ab1dc0SHong Zhang       ierr = PetscLLCondensedAddSorted(anz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1343f96d37fcSHong Zhang     }
13444e938277SHong Zhang 
1345f96d37fcSHong Zhang     /* add received col data into lnk */
1346f96d37fcSHong Zhang     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1347f96d37fcSHong Zhang       if (i == *nextrow[k]) { /* i-th row */
1348f96d37fcSHong Zhang         nzi  = *(nextci[k]+1) - *nextci[k];
1349f96d37fcSHong Zhang         Jptr = buf_rj[k] + *nextci[k];
13504e938277SHong Zhang         ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1351f96d37fcSHong Zhang         nextrow[k]++; nextci[k]++;
1352f96d37fcSHong Zhang       }
1353f96d37fcSHong Zhang     }
13544e938277SHong Zhang     nnz = lnk[0];
1355f96d37fcSHong Zhang 
1356f96d37fcSHong Zhang     /* if free space is not available, make more free space */
1357f96d37fcSHong Zhang     if (current_space->local_remaining<nnz) {
1358f91af8c7SBarry Smith       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(nnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
1359f96d37fcSHong Zhang       nspacedouble++;
1360f96d37fcSHong Zhang     }
1361f96d37fcSHong Zhang     /* copy data into free space, then initialize lnk */
13624e938277SHong Zhang     ierr = PetscLLCondensedClean(aN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
1363f96d37fcSHong Zhang     ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr);
13642205254eSKarl Rupp 
1365f96d37fcSHong Zhang     current_space->array           += nnz;
1366f96d37fcSHong Zhang     current_space->local_used      += nnz;
1367f96d37fcSHong Zhang     current_space->local_remaining -= nnz;
13682205254eSKarl Rupp 
1369f96d37fcSHong Zhang     bi[i+1] = bi[i] + nnz;
1370f96d37fcSHong Zhang     if (nnz > rmax) rmax = nnz;
1371f96d37fcSHong Zhang   }
1372f96d37fcSHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1373f96d37fcSHong Zhang 
1374854ce69bSBarry Smith   ierr = PetscMalloc1(bi[pn]+1,&bj);CHKERRQ(ierr);
1375f96d37fcSHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
13762205254eSKarl Rupp 
1377118727c9SMark F. Adams   afill_tmp = (PetscReal)bi[pn]/(pdti[pn] + poti[pon] + ai[am]+1);
1378f96d37fcSHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
1379d6ab1dc0SHong Zhang   ierr = PetscLLCondensedDestroy(lnk,lnkbt);CHKERRQ(ierr);
13804e938277SHong Zhang   ierr = MatDestroy(&POt);CHKERRQ(ierr);
13814e938277SHong Zhang   ierr = MatDestroy(&PDt);CHKERRQ(ierr);
1382f96d37fcSHong Zhang 
13831c7d5954SHong Zhang   /* create symbolic parallel matrix Cmpi - why cannot be assembled in Numeric part   */
13841c7d5954SHong Zhang   /*----------------------------------------------------------------------------------*/
13851795a4d1SJed Brown   ierr = PetscCalloc1(rmax+1,&vals);CHKERRQ(ierr);
1386f96d37fcSHong Zhang 
1387f96d37fcSHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
1388f96d37fcSHong Zhang   ierr = MatSetSizes(Cmpi,pn,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
138933d57670SJed Brown   ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(A->cmap->bs));CHKERRQ(ierr);
1390f96d37fcSHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
1391f96d37fcSHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
1392f96d37fcSHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
1393f96d37fcSHong Zhang   ierr = MatSetBlockSize(Cmpi,1);CHKERRQ(ierr);
1394f96d37fcSHong Zhang   for (i=0; i<pn; i++) {
1395f96d37fcSHong Zhang     row  = i + rstart;
1396f96d37fcSHong Zhang     nnz  = bi[i+1] - bi[i];
1397f96d37fcSHong Zhang     Jptr = bj + bi[i];
1398f96d37fcSHong Zhang     ierr = MatSetValues(Cmpi,1,&row,nnz,Jptr,vals,INSERT_VALUES);CHKERRQ(ierr);
1399f96d37fcSHong Zhang   }
1400f96d37fcSHong Zhang   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1401f96d37fcSHong Zhang   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1402f96d37fcSHong Zhang   ierr = PetscFree(vals);CHKERRQ(ierr);
1403f96d37fcSHong Zhang 
1404f96d37fcSHong Zhang   merge->bi        = bi;
1405f96d37fcSHong Zhang   merge->bj        = bj;
1406f96d37fcSHong Zhang   merge->coi       = coi;
1407f96d37fcSHong Zhang   merge->coj       = coj;
1408f96d37fcSHong Zhang   merge->buf_ri    = buf_ri;
1409f96d37fcSHong Zhang   merge->buf_rj    = buf_rj;
1410f96d37fcSHong Zhang   merge->owners_co = owners_co;
1411f96d37fcSHong Zhang 
1412f96d37fcSHong Zhang   /* attach the supporting struct to Cmpi for reuse */
1413f96d37fcSHong Zhang   c           = (Mat_MPIAIJ*)Cmpi->data;
1414f96d37fcSHong Zhang   c->ptap     = ptap;
14150298fd71SBarry Smith   ptap->api   = NULL;
14160298fd71SBarry Smith   ptap->apj   = NULL;
1417f96d37fcSHong Zhang   ptap->merge = merge;
1418375ed354SHong Zhang   ptap->destroy   = Cmpi->ops->destroy;
1419375ed354SHong Zhang   ptap->duplicate = Cmpi->ops->duplicate;
1420375ed354SHong Zhang 
1421375ed354SHong Zhang   Cmpi->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable;
1422375ed354SHong Zhang   Cmpi->ops->destroy                 = MatDestroy_MPIAIJ_PtAP;
1423375ed354SHong Zhang   Cmpi->ops->duplicate               = MatDuplicate_MPIAIJ_MatPtAP;
1424d6ab1dc0SHong Zhang 
1425d6ab1dc0SHong Zhang   *C = Cmpi;
1426d6ab1dc0SHong Zhang #if defined(PETSC_USE_INFO)
1427d6ab1dc0SHong Zhang   if (bi[pn] != 0) {
142857622a8eSBarry Smith     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
142957622a8eSBarry Smith     ierr = PetscInfo1(Cmpi,"Use MatTransposeMatMult(A,B,MatReuse,%g,&C) for best performance.\n",(double)afill);CHKERRQ(ierr);
1430d6ab1dc0SHong Zhang   } else {
1431d6ab1dc0SHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
1432d6ab1dc0SHong Zhang   }
1433d6ab1dc0SHong Zhang #endif
1434d6ab1dc0SHong Zhang   PetscFunctionReturn(0);
1435d6ab1dc0SHong Zhang }
1436d6ab1dc0SHong Zhang 
1437d6ab1dc0SHong Zhang #undef __FUNCT__
14386da69ca6SHong Zhang #define __FUNCT__ "MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ"
14396da69ca6SHong Zhang PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ(Mat P,Mat A,Mat C)
1440d6ab1dc0SHong Zhang {
1441d6ab1dc0SHong Zhang   PetscErrorCode      ierr;
1442d6ab1dc0SHong Zhang   Mat_Merge_SeqsToMPI *merge;
1443d6ab1dc0SHong Zhang   Mat_MPIAIJ          *p =(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
1444d6ab1dc0SHong Zhang   Mat_SeqAIJ          *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data;
1445d6ab1dc0SHong Zhang   Mat_PtAPMPI         *ptap;
1446e745005bSHong Zhang   PetscInt            *adj;
1447e745005bSHong Zhang   PetscInt            i,j,k,anz,pnz,row,*cj,nexta;
1448e745005bSHong Zhang   MatScalar           *ada,*ca,valtmp;
1449d6ab1dc0SHong Zhang   PetscInt            am  =A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n;
1450ce94432eSBarry Smith   MPI_Comm            comm;
1451d6ab1dc0SHong Zhang   PetscMPIInt         size,rank,taga,*len_s;
1452d6ab1dc0SHong Zhang   PetscInt            *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci;
1453d6ab1dc0SHong Zhang   PetscInt            **buf_ri,**buf_rj;
1454d6ab1dc0SHong Zhang   PetscInt            cnz=0,*bj_i,*bi,*bj,bnz,nextcj;  /* bi,bj,ba: local array of C(mpi mat) */
1455d6ab1dc0SHong Zhang   MPI_Request         *s_waits,*r_waits;
1456d6ab1dc0SHong Zhang   MPI_Status          *status;
1457d6ab1dc0SHong Zhang   MatScalar           **abuf_r,*ba_i,*pA,*coa,*ba;
1458d6ab1dc0SHong Zhang   PetscInt            *ai,*aj,*coi,*coj;
1459a2ea699eSBarry Smith   PetscInt            *poJ,*pdJ;
1460d6ab1dc0SHong Zhang   Mat                 A_loc;
1461d6ab1dc0SHong Zhang   Mat_SeqAIJ          *a_loc;
1462d6ab1dc0SHong Zhang 
1463d6ab1dc0SHong Zhang   PetscFunctionBegin;
1464ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr);
1465d6ab1dc0SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1466d6ab1dc0SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1467d6ab1dc0SHong Zhang 
1468d6ab1dc0SHong Zhang   ptap  = c->ptap;
1469d6ab1dc0SHong Zhang   merge = ptap->merge;
1470d6ab1dc0SHong Zhang 
1471e745005bSHong Zhang   /* 2) compute numeric C_seq = P_loc^T*A_loc */
1472e745005bSHong Zhang   /*------------------------------------------*/
1473d6ab1dc0SHong Zhang   /* get data from symbolic products */
1474d6ab1dc0SHong Zhang   coi    = merge->coi; coj = merge->coj;
1475854ce69bSBarry Smith   ierr   = PetscCalloc1(coi[pon]+1,&coa);CHKERRQ(ierr);
1476d6ab1dc0SHong Zhang   bi     = merge->bi; bj = merge->bj;
1477d6ab1dc0SHong Zhang   owners = merge->rowmap->range;
14781795a4d1SJed Brown   ierr   = PetscCalloc1(bi[cm]+1,&ba);CHKERRQ(ierr);
1479d6ab1dc0SHong Zhang 
1480d6ab1dc0SHong Zhang   /* get A_loc by taking all local rows of A */
1481d6ab1dc0SHong Zhang   A_loc = ptap->A_loc;
1482d6ab1dc0SHong Zhang   ierr  = MatMPIAIJGetLocalMat(A,MAT_REUSE_MATRIX,&A_loc);CHKERRQ(ierr);
1483d6ab1dc0SHong Zhang   a_loc = (Mat_SeqAIJ*)(A_loc)->data;
1484d6ab1dc0SHong Zhang   ai    = a_loc->i;
1485d6ab1dc0SHong Zhang   aj    = a_loc->j;
1486d6ab1dc0SHong Zhang 
1487d6ab1dc0SHong Zhang   for (i=0; i<am; i++) {
1488d6ab1dc0SHong Zhang     anz = ai[i+1] - ai[i];
1489d6ab1dc0SHong Zhang     adj = aj + ai[i];
1490d6ab1dc0SHong Zhang     ada = a_loc->a + ai[i];
1491d6ab1dc0SHong Zhang 
1492d6ab1dc0SHong Zhang     /* 2-b) Compute Cseq = P_loc[i,:]^T*A[i,:] using outer product */
1493e745005bSHong Zhang     /*-------------------------------------------------------------*/
1494d6ab1dc0SHong Zhang     /* put the value into Co=(p->B)^T*A (off-diagonal part, send to others) */
1495d6ab1dc0SHong Zhang     pnz = po->i[i+1] - po->i[i];
1496d6ab1dc0SHong Zhang     poJ = po->j + po->i[i];
1497d6ab1dc0SHong Zhang     pA  = po->a + po->i[i];
1498d6ab1dc0SHong Zhang     for (j=0; j<pnz; j++) {
1499d6ab1dc0SHong Zhang       row = poJ[j];
1500d6ab1dc0SHong Zhang       cj  = coj + coi[row];
1501d6ab1dc0SHong Zhang       ca  = coa + coi[row];
1502e745005bSHong Zhang       /* perform sparse axpy */
1503e745005bSHong Zhang       nexta  = 0;
1504d6ab1dc0SHong Zhang       valtmp = pA[j];
1505e745005bSHong Zhang       for (k=0; nexta<anz; k++) {
1506e745005bSHong Zhang         if (cj[k] == adj[nexta]) {
1507e745005bSHong Zhang           ca[k] += valtmp*ada[nexta];
1508e745005bSHong Zhang           nexta++;
1509d6ab1dc0SHong Zhang         }
1510e745005bSHong Zhang       }
1511e745005bSHong Zhang       ierr = PetscLogFlops(2.0*anz);CHKERRQ(ierr);
1512d6ab1dc0SHong Zhang     }
1513d6ab1dc0SHong Zhang 
1514d6ab1dc0SHong Zhang     /* put the value into Cd (diagonal part) */
1515d6ab1dc0SHong Zhang     pnz = pd->i[i+1] - pd->i[i];
1516d6ab1dc0SHong Zhang     pdJ = pd->j + pd->i[i];
1517d6ab1dc0SHong Zhang     pA  = pd->a + pd->i[i];
1518d6ab1dc0SHong Zhang     for (j=0; j<pnz; j++) {
1519d6ab1dc0SHong Zhang       row = pdJ[j];
1520d6ab1dc0SHong Zhang       cj  = bj + bi[row];
1521d6ab1dc0SHong Zhang       ca  = ba + bi[row];
1522e745005bSHong Zhang       /* perform sparse axpy */
1523e745005bSHong Zhang       nexta  = 0;
1524d6ab1dc0SHong Zhang       valtmp = pA[j];
1525e745005bSHong Zhang       for (k=0; nexta<anz; k++) {
1526e745005bSHong Zhang         if (cj[k] == adj[nexta]) {
1527e745005bSHong Zhang           ca[k] += valtmp*ada[nexta];
1528e745005bSHong Zhang           nexta++;
1529d6ab1dc0SHong Zhang         }
1530e745005bSHong Zhang       }
1531e745005bSHong Zhang       ierr = PetscLogFlops(2.0*anz);CHKERRQ(ierr);
1532d6ab1dc0SHong Zhang     }
1533d6ab1dc0SHong Zhang   }
1534d6ab1dc0SHong Zhang 
1535d6ab1dc0SHong Zhang   /* 3) send and recv matrix values coa */
1536d6ab1dc0SHong Zhang   /*------------------------------------*/
1537d6ab1dc0SHong Zhang   buf_ri = merge->buf_ri;
1538d6ab1dc0SHong Zhang   buf_rj = merge->buf_rj;
1539d6ab1dc0SHong Zhang   len_s  = merge->len_s;
1540d6ab1dc0SHong Zhang   ierr   = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr);
1541d6ab1dc0SHong Zhang   ierr   = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr);
1542d6ab1dc0SHong Zhang 
1543dcca6d9dSJed Brown   ierr = PetscMalloc2(merge->nsend+1,&s_waits,size,&status);CHKERRQ(ierr);
1544d6ab1dc0SHong Zhang   for (proc=0,k=0; proc<size; proc++) {
1545d6ab1dc0SHong Zhang     if (!len_s[proc]) continue;
1546d6ab1dc0SHong Zhang     i    = merge->owners_co[proc];
1547d6ab1dc0SHong Zhang     ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr);
1548d6ab1dc0SHong Zhang     k++;
1549d6ab1dc0SHong Zhang   }
1550d6ab1dc0SHong Zhang   if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);}
1551d6ab1dc0SHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);}
1552d6ab1dc0SHong Zhang 
1553d6ab1dc0SHong Zhang   ierr = PetscFree2(s_waits,status);CHKERRQ(ierr);
1554d6ab1dc0SHong Zhang   ierr = PetscFree(r_waits);CHKERRQ(ierr);
1555d6ab1dc0SHong Zhang   ierr = PetscFree(coa);CHKERRQ(ierr);
1556d6ab1dc0SHong Zhang 
1557d6ab1dc0SHong Zhang   /* 4) insert local Cseq and received values into Cmpi */
1558d6ab1dc0SHong Zhang   /*----------------------------------------------------*/
1559dcca6d9dSJed Brown   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
1560d6ab1dc0SHong Zhang   for (k=0; k<merge->nrecv; k++) {
1561e745005bSHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1562d6ab1dc0SHong Zhang     nrows       = *(buf_ri_k[k]);
1563d6ab1dc0SHong Zhang     nextrow[k]  = buf_ri_k[k]+1;  /* next row number of k-th recved i-structure */
1564d6ab1dc0SHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
1565d6ab1dc0SHong Zhang   }
1566d6ab1dc0SHong Zhang 
1567d6ab1dc0SHong Zhang   for (i=0; i<cm; i++) {
1568d6ab1dc0SHong Zhang     row  = owners[rank] + i; /* global row index of C_seq */
1569d6ab1dc0SHong Zhang     bj_i = bj + bi[i];  /* col indices of the i-th row of C */
1570d6ab1dc0SHong Zhang     ba_i = ba + bi[i];
1571d6ab1dc0SHong Zhang     bnz  = bi[i+1] - bi[i];
1572d6ab1dc0SHong Zhang     /* add received vals into ba */
1573d6ab1dc0SHong Zhang     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1574d6ab1dc0SHong Zhang       /* i-th row */
1575d6ab1dc0SHong Zhang       if (i == *nextrow[k]) {
1576d6ab1dc0SHong Zhang         cnz    = *(nextci[k]+1) - *nextci[k];
1577d6ab1dc0SHong Zhang         cj     = buf_rj[k] + *(nextci[k]);
1578d6ab1dc0SHong Zhang         ca     = abuf_r[k] + *(nextci[k]);
1579d6ab1dc0SHong Zhang         nextcj = 0;
1580d6ab1dc0SHong Zhang         for (j=0; nextcj<cnz; j++) {
1581d6ab1dc0SHong Zhang           if (bj_i[j] == cj[nextcj]) { /* bcol == ccol */
1582d6ab1dc0SHong Zhang             ba_i[j] += ca[nextcj++];
1583d6ab1dc0SHong Zhang           }
1584d6ab1dc0SHong Zhang         }
1585d6ab1dc0SHong Zhang         nextrow[k]++; nextci[k]++;
1586d6ab1dc0SHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1587d6ab1dc0SHong Zhang       }
1588d6ab1dc0SHong Zhang     }
1589d6ab1dc0SHong Zhang     ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr);
1590d6ab1dc0SHong Zhang   }
1591d6ab1dc0SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1592d6ab1dc0SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1593d6ab1dc0SHong Zhang 
1594d6ab1dc0SHong Zhang   ierr = PetscFree(ba);CHKERRQ(ierr);
1595d6ab1dc0SHong Zhang   ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr);
1596d6ab1dc0SHong Zhang   ierr = PetscFree(abuf_r);CHKERRQ(ierr);
1597d6ab1dc0SHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1598d6ab1dc0SHong Zhang   PetscFunctionReturn(0);
1599d6ab1dc0SHong Zhang }
1600d6ab1dc0SHong Zhang 
1601c9ba551fSBarry Smith PetscErrorCode MatDuplicate_MPIAIJ_MatPtAP(Mat, MatDuplicateOption,Mat*);
16022bbb1c24SHong Zhang /* This routine is modified from MatPtAPSymbolic_MPIAIJ_MPIAIJ();
16032bbb1c24SHong Zhang    differ from MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable in using LLCondensedCreate_Scalable() */
1604d6ab1dc0SHong Zhang #undef __FUNCT__
16056da69ca6SHong Zhang #define __FUNCT__ "MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ"
16066da69ca6SHong Zhang PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(Mat P,Mat A,PetscReal fill,Mat *C)
1607d6ab1dc0SHong Zhang {
1608d6ab1dc0SHong Zhang   PetscErrorCode      ierr;
1609d6ab1dc0SHong Zhang   Mat                 Cmpi,A_loc,POt,PDt;
1610d6ab1dc0SHong Zhang   Mat_PtAPMPI         *ptap;
16110298fd71SBarry Smith   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
1612*33ba5e2fSHong Zhang   Mat_MPIAIJ          *p=(Mat_MPIAIJ*)P->data,*a=(Mat_MPIAIJ*)A->data,*c;
1613d6ab1dc0SHong Zhang   PetscInt            *pdti,*pdtj,*poti,*potj,*ptJ;
1614d6ab1dc0SHong Zhang   PetscInt            nnz;
1615d6ab1dc0SHong Zhang   PetscInt            *lnk,*owners_co,*coi,*coj,i,k,pnz,row;
1616d6ab1dc0SHong Zhang   PetscInt            am  =A->rmap->n,pn=P->cmap->n;
1617ce94432eSBarry Smith   MPI_Comm            comm;
1618d6ab1dc0SHong Zhang   PetscMPIInt         size,rank,tagi,tagj,*len_si,*len_s,*len_ri;
1619d6ab1dc0SHong Zhang   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
1620d6ab1dc0SHong Zhang   PetscInt            len,proc,*dnz,*onz,*owners;
1621d6ab1dc0SHong Zhang   PetscInt            nzi,*bi,*bj;
1622d6ab1dc0SHong Zhang   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
1623d6ab1dc0SHong Zhang   MPI_Request         *swaits,*rwaits;
1624d6ab1dc0SHong Zhang   MPI_Status          *sstatus,rstatus;
1625d6ab1dc0SHong Zhang   Mat_Merge_SeqsToMPI *merge;
1626d6ab1dc0SHong Zhang   PetscInt            *ai,*aj,*Jptr,anz,*prmap=p->garray,pon,nspacedouble=0,j;
1627d6ab1dc0SHong Zhang   PetscReal           afill  =1.0,afill_tmp;
16284b38b95cSHong Zhang   PetscInt            rstart = P->cmap->rstart,rmax,aN=A->cmap->N,Armax;
1629d6ab1dc0SHong Zhang   PetscScalar         *vals;
1630d6ab1dc0SHong Zhang   Mat_SeqAIJ          *a_loc,*pdt,*pot;
16310acc65b4SHong Zhang   PetscTable          ta;
1632d6ab1dc0SHong Zhang 
1633d6ab1dc0SHong Zhang   PetscFunctionBegin;
1634ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
1635d6ab1dc0SHong Zhang   /* check if matrix local sizes are compatible */
16366c4ed002SBarry Smith   if (A->rmap->rstart != P->rmap->rstart || A->rmap->rend != P->rmap->rend) SETERRQ4(comm,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, A (%D, %D) != P (%D,%D)",A->rmap->rstart,A->rmap->rend,P->rmap->rstart,P->rmap->rend);
1637d6ab1dc0SHong Zhang 
1638d6ab1dc0SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1639d6ab1dc0SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1640d6ab1dc0SHong Zhang 
1641d6ab1dc0SHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
1642b00a9115SJed Brown   ierr = PetscNew(&ptap);CHKERRQ(ierr);
1643d6ab1dc0SHong Zhang 
1644d6ab1dc0SHong Zhang   /* get A_loc by taking all local rows of A */
1645d6ab1dc0SHong Zhang   ierr = MatMPIAIJGetLocalMat(A,MAT_INITIAL_MATRIX,&A_loc);CHKERRQ(ierr);
16462205254eSKarl Rupp 
1647d6ab1dc0SHong Zhang   ptap->A_loc = A_loc;
1648d6ab1dc0SHong Zhang   a_loc       = (Mat_SeqAIJ*)(A_loc)->data;
1649d6ab1dc0SHong Zhang   ai          = a_loc->i;
1650d6ab1dc0SHong Zhang   aj          = a_loc->j;
1651d6ab1dc0SHong Zhang 
1652d6ab1dc0SHong Zhang   /* determine symbolic Co=(p->B)^T*A - send to others */
1653d6ab1dc0SHong Zhang   /*----------------------------------------------------*/
1654d6ab1dc0SHong Zhang   ierr = MatTransposeSymbolic_SeqAIJ(p->A,&PDt);CHKERRQ(ierr);
1655d6ab1dc0SHong Zhang   pdt  = (Mat_SeqAIJ*)PDt->data;
1656d6ab1dc0SHong Zhang   pdti = pdt->i; pdtj = pdt->j;
1657d6ab1dc0SHong Zhang 
1658d6ab1dc0SHong Zhang   ierr = MatTransposeSymbolic_SeqAIJ(p->B,&POt);CHKERRQ(ierr);
1659d6ab1dc0SHong Zhang   pot  = (Mat_SeqAIJ*)POt->data;
1660d6ab1dc0SHong Zhang   poti = pot->i; potj = pot->j;
1661d6ab1dc0SHong Zhang 
1662d6ab1dc0SHong Zhang   /* then, compute symbolic Co = (p->B)^T*A */
1663d6ab1dc0SHong Zhang   pon    = (p->B)->cmap->n; /* total num of rows to be sent to other processors
1664d6ab1dc0SHong Zhang                          >= (num of nonzero rows of C_seq) - pn */
1665854ce69bSBarry Smith   ierr   = PetscMalloc1(pon+1,&coi);CHKERRQ(ierr);
1666d6ab1dc0SHong Zhang   coi[0] = 0;
1667d6ab1dc0SHong Zhang 
1668d6ab1dc0SHong Zhang   /* set initial free space to be fill*(nnz(p->B) + nnz(A)) */
1669f91af8c7SBarry Smith   nnz           = PetscRealIntMultTruncate(fill,PetscIntSumTruncate(poti[pon],ai[am]));
1670a2ea699eSBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
1671d6ab1dc0SHong Zhang   current_space = free_space;
167219f0e57aSHong Zhang 
1673d6ab1dc0SHong Zhang   /* create and initialize a linked list */
1674*33ba5e2fSHong Zhang   ierr = PetscTableCreate(A->cmap->n + a->B->cmap->N,aN,&ta);CHKERRQ(ierr);
16754b38b95cSHong Zhang   MatRowMergeMax_SeqAIJ(a_loc,am,ta);
16760acc65b4SHong Zhang   ierr = PetscTableGetCount(ta,&Armax);CHKERRQ(ierr);
1677*33ba5e2fSHong Zhang 
16780acc65b4SHong Zhang   ierr = PetscLLCondensedCreate_Scalable(Armax,&lnk);CHKERRQ(ierr);
1679d6ab1dc0SHong Zhang 
1680d6ab1dc0SHong Zhang   for (i=0; i<pon; i++) {
1681d6ab1dc0SHong Zhang     pnz = poti[i+1] - poti[i];
1682d6ab1dc0SHong Zhang     ptJ = potj + poti[i];
1683d6ab1dc0SHong Zhang     for (j=0; j<pnz; j++) {
1684d6ab1dc0SHong Zhang       row  = ptJ[j]; /* row of A_loc == col of Pot */
1685d6ab1dc0SHong Zhang       anz  = ai[row+1] - ai[row];
1686d6ab1dc0SHong Zhang       Jptr = aj + ai[row];
1687d6ab1dc0SHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
1688d6ab1dc0SHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(anz,Jptr,lnk);CHKERRQ(ierr);
1689d6ab1dc0SHong Zhang     }
1690d6ab1dc0SHong Zhang     nnz = lnk[0];
1691d6ab1dc0SHong Zhang 
1692d6ab1dc0SHong Zhang     /* If free space is not available, double the total space in the list */
1693d6ab1dc0SHong Zhang     if (current_space->local_remaining<nnz) {
1694f91af8c7SBarry Smith       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(nnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
1695d6ab1dc0SHong Zhang       nspacedouble++;
1696d6ab1dc0SHong Zhang     }
1697d6ab1dc0SHong Zhang 
1698d6ab1dc0SHong Zhang     /* Copy data into free space, and zero out denserows */
1699d6ab1dc0SHong Zhang     ierr = PetscLLCondensedClean_Scalable(nnz,current_space->array,lnk);CHKERRQ(ierr);
17002205254eSKarl Rupp 
1701d6ab1dc0SHong Zhang     current_space->array           += nnz;
1702d6ab1dc0SHong Zhang     current_space->local_used      += nnz;
1703d6ab1dc0SHong Zhang     current_space->local_remaining -= nnz;
17042205254eSKarl Rupp 
1705d6ab1dc0SHong Zhang     coi[i+1] = coi[i] + nnz;
1706d6ab1dc0SHong Zhang   }
1707d6ab1dc0SHong Zhang 
1708854ce69bSBarry Smith   ierr = PetscMalloc1(coi[pon]+1,&coj);CHKERRQ(ierr);
1709d6ab1dc0SHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr);
17100acc65b4SHong Zhang   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr); /* must destroy to get a new one for C */
17112205254eSKarl Rupp 
1712118727c9SMark F. Adams   afill_tmp = (PetscReal)coi[pon]/(poti[pon] + ai[am]+1);
1713d6ab1dc0SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
1714d6ab1dc0SHong Zhang 
1715d6ab1dc0SHong Zhang   /* send j-array (coj) of Co to other processors */
1716d6ab1dc0SHong Zhang   /*----------------------------------------------*/
1717d6ab1dc0SHong Zhang   /* determine row ownership */
1718b00a9115SJed Brown   ierr = PetscNew(&merge);CHKERRQ(ierr);
1719d6ab1dc0SHong Zhang   ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr);
17202205254eSKarl Rupp 
1721d6ab1dc0SHong Zhang   merge->rowmap->n  = pn;
1722d6ab1dc0SHong Zhang   merge->rowmap->bs = 1;
17232205254eSKarl Rupp 
1724d6ab1dc0SHong Zhang   ierr   = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr);
1725d6ab1dc0SHong Zhang   owners = merge->rowmap->range;
1726d6ab1dc0SHong Zhang 
1727d6ab1dc0SHong Zhang   /* determine the number of messages to send, their lengths */
17281795a4d1SJed Brown   ierr = PetscCalloc1(size,&len_si);CHKERRQ(ierr);
1729785e854fSJed Brown   ierr = PetscMalloc1(size,&merge->len_s);CHKERRQ(ierr);
17302205254eSKarl Rupp 
1731d6ab1dc0SHong Zhang   len_s        = merge->len_s;
1732d6ab1dc0SHong Zhang   merge->nsend = 0;
1733d6ab1dc0SHong Zhang 
1734854ce69bSBarry Smith   ierr = PetscMalloc1(size+2,&owners_co);CHKERRQ(ierr);
1735d6ab1dc0SHong Zhang   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
1736d6ab1dc0SHong Zhang 
1737d6ab1dc0SHong Zhang   proc = 0;
1738d6ab1dc0SHong Zhang   for (i=0; i<pon; i++) {
1739d6ab1dc0SHong Zhang     while (prmap[i] >= owners[proc+1]) proc++;
1740d6ab1dc0SHong Zhang     len_si[proc]++;  /* num of rows in Co to be sent to [proc] */
1741d6ab1dc0SHong Zhang     len_s[proc] += coi[i+1] - coi[i];
1742d6ab1dc0SHong Zhang   }
1743d6ab1dc0SHong Zhang 
1744d6ab1dc0SHong Zhang   len          = 0; /* max length of buf_si[] */
1745d6ab1dc0SHong Zhang   owners_co[0] = 0;
1746d6ab1dc0SHong Zhang   for (proc=0; proc<size; proc++) {
1747d6ab1dc0SHong Zhang     owners_co[proc+1] = owners_co[proc] + len_si[proc];
1748d6ab1dc0SHong Zhang     if (len_si[proc]) {
1749d6ab1dc0SHong Zhang       merge->nsend++;
1750d6ab1dc0SHong Zhang       len_si[proc] = 2*(len_si[proc] + 1);
1751d6ab1dc0SHong Zhang       len         += len_si[proc];
1752d6ab1dc0SHong Zhang     }
1753d6ab1dc0SHong Zhang   }
1754d6ab1dc0SHong Zhang 
1755d6ab1dc0SHong Zhang   /* determine the number and length of messages to receive for coi and coj  */
17560298fd71SBarry Smith   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&merge->nrecv);CHKERRQ(ierr);
1757d6ab1dc0SHong Zhang   ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr);
1758d6ab1dc0SHong Zhang 
1759d6ab1dc0SHong Zhang   /* post the Irecv and Isend of coj */
1760d6ab1dc0SHong Zhang   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
1761d6ab1dc0SHong Zhang   ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
1762854ce69bSBarry Smith   ierr = PetscMalloc1(merge->nsend+1,&swaits);CHKERRQ(ierr);
1763d6ab1dc0SHong Zhang   for (proc=0, k=0; proc<size; proc++) {
1764d6ab1dc0SHong Zhang     if (!len_s[proc]) continue;
1765d6ab1dc0SHong Zhang     i    = owners_co[proc];
1766d6ab1dc0SHong Zhang     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
1767d6ab1dc0SHong Zhang     k++;
1768d6ab1dc0SHong Zhang   }
1769d6ab1dc0SHong Zhang 
1770d6ab1dc0SHong Zhang   /* receives and sends of coj are complete */
1771785e854fSJed Brown   ierr = PetscMalloc1(size,&sstatus);CHKERRQ(ierr);
1772d6ab1dc0SHong Zhang   for (i=0; i<merge->nrecv; i++) {
1773c280ed6aSJed Brown     PetscMPIInt icompleted;
1774c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1775d6ab1dc0SHong Zhang   }
1776d6ab1dc0SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1777d6ab1dc0SHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1778d6ab1dc0SHong Zhang 
17790acc65b4SHong Zhang   /* add received column indices into table to update Armax */
1780*33ba5e2fSHong Zhang   /* Armax can be as large as aN if a P[row,:] is dense, see src/ksp/ksp/examples/tutorials/ex56.c! */
17810acc65b4SHong Zhang   for (k=0; k<merge->nrecv; k++) {/* k-th received message */
17820acc65b4SHong Zhang     Jptr = buf_rj[k];
17830acc65b4SHong Zhang     for (j=0; j<merge->len_r[k]; j++) {
17840acc65b4SHong Zhang       ierr = PetscTableAdd(ta,*(Jptr+j)+1,1,INSERT_VALUES);CHKERRQ(ierr);
17850acc65b4SHong Zhang     }
17860acc65b4SHong Zhang   }
17870acc65b4SHong Zhang   ierr = PetscTableGetCount(ta,&Armax);CHKERRQ(ierr);
1788*33ba5e2fSHong Zhang   /* printf("Armax %d, an %d + Bn %d = %d, aN %d\n",Armax,A->cmap->n,a->B->cmap->N,A->cmap->n+a->B->cmap->N,aN); */
17890acc65b4SHong Zhang 
1790d6ab1dc0SHong Zhang   /* send and recv coi */
1791d6ab1dc0SHong Zhang   /*-------------------*/
1792d6ab1dc0SHong Zhang   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
1793d6ab1dc0SHong Zhang   ierr   = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
1794854ce69bSBarry Smith   ierr   = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr);
1795d6ab1dc0SHong Zhang   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
1796d6ab1dc0SHong Zhang   for (proc=0,k=0; proc<size; proc++) {
1797d6ab1dc0SHong Zhang     if (!len_s[proc]) continue;
1798d6ab1dc0SHong Zhang     /* form outgoing message for i-structure:
1799d6ab1dc0SHong Zhang          buf_si[0]:                 nrows to be sent
1800d6ab1dc0SHong Zhang                [1:nrows]:           row index (global)
1801d6ab1dc0SHong Zhang                [nrows+1:2*nrows+1]: i-structure index
1802d6ab1dc0SHong Zhang     */
1803d6ab1dc0SHong Zhang     /*-------------------------------------------*/
1804d6ab1dc0SHong Zhang     nrows       = len_si[proc]/2 - 1;
1805d6ab1dc0SHong Zhang     buf_si_i    = buf_si + nrows+1;
1806d6ab1dc0SHong Zhang     buf_si[0]   = nrows;
1807d6ab1dc0SHong Zhang     buf_si_i[0] = 0;
1808d6ab1dc0SHong Zhang     nrows       = 0;
1809d6ab1dc0SHong Zhang     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
1810d6ab1dc0SHong Zhang       nzi               = coi[i+1] - coi[i];
1811d6ab1dc0SHong Zhang       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi;  /* i-structure */
1812d6ab1dc0SHong Zhang       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
1813d6ab1dc0SHong Zhang       nrows++;
1814d6ab1dc0SHong Zhang     }
1815d6ab1dc0SHong Zhang     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
1816d6ab1dc0SHong Zhang     k++;
1817d6ab1dc0SHong Zhang     buf_si += len_si[proc];
1818d6ab1dc0SHong Zhang   }
1819d6ab1dc0SHong Zhang   i = merge->nrecv;
1820d6ab1dc0SHong Zhang   while (i--) {
1821c280ed6aSJed Brown     PetscMPIInt icompleted;
1822c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1823d6ab1dc0SHong Zhang   }
1824d6ab1dc0SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1825d6ab1dc0SHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1826d6ab1dc0SHong Zhang   ierr = PetscFree(len_si);CHKERRQ(ierr);
1827d6ab1dc0SHong Zhang   ierr = PetscFree(len_ri);CHKERRQ(ierr);
1828d6ab1dc0SHong Zhang   ierr = PetscFree(swaits);CHKERRQ(ierr);
1829d6ab1dc0SHong Zhang   ierr = PetscFree(sstatus);CHKERRQ(ierr);
1830d6ab1dc0SHong Zhang   ierr = PetscFree(buf_s);CHKERRQ(ierr);
1831d6ab1dc0SHong Zhang 
1832d6ab1dc0SHong Zhang   /* compute the local portion of C (mpi mat) */
1833d6ab1dc0SHong Zhang   /*------------------------------------------*/
1834d6ab1dc0SHong Zhang   /* allocate bi array and free space for accumulating nonzero column info */
1835854ce69bSBarry Smith   ierr  = PetscMalloc1(pn+1,&bi);CHKERRQ(ierr);
1836d6ab1dc0SHong Zhang   bi[0] = 0;
1837d6ab1dc0SHong Zhang 
1838d6ab1dc0SHong Zhang   /* set initial free space to be fill*(nnz(P) + nnz(AP)) */
1839f91af8c7SBarry Smith   nnz           = PetscRealIntMultTruncate(fill,PetscIntSumTruncate(pdti[pn],PetscIntSumTruncate(poti[pon],ai[am])));
1840a2ea699eSBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
1841d6ab1dc0SHong Zhang   current_space = free_space;
1842d6ab1dc0SHong Zhang 
1843dcca6d9dSJed Brown   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
1844d6ab1dc0SHong Zhang   for (k=0; k<merge->nrecv; k++) {
1845d6ab1dc0SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1846d6ab1dc0SHong Zhang     nrows       = *buf_ri_k[k];
1847d6ab1dc0SHong Zhang     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
18482205254eSKarl Rupp     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* points to the next i-structure of k-th recieved i-structure  */
1849d6ab1dc0SHong Zhang   }
1850d6ab1dc0SHong Zhang 
18510acc65b4SHong Zhang   ierr = PetscLLCondensedCreate_Scalable(Armax,&lnk);CHKERRQ(ierr);
1852d6ab1dc0SHong Zhang   ierr = MatPreallocateInitialize(comm,pn,A->cmap->n,dnz,onz);CHKERRQ(ierr);
1853d6ab1dc0SHong Zhang   rmax = 0;
1854d6ab1dc0SHong Zhang   for (i=0; i<pn; i++) {
1855d6ab1dc0SHong Zhang     /* add pdt[i,:]*AP into lnk */
1856d6ab1dc0SHong Zhang     pnz = pdti[i+1] - pdti[i];
1857d6ab1dc0SHong Zhang     ptJ = pdtj + pdti[i];
1858d6ab1dc0SHong Zhang     for (j=0; j<pnz; j++) {
1859d6ab1dc0SHong Zhang       row  = ptJ[j];  /* row of AP == col of Pt */
1860d6ab1dc0SHong Zhang       anz  = ai[row+1] - ai[row];
1861d6ab1dc0SHong Zhang       Jptr = aj + ai[row];
1862d6ab1dc0SHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
1863d6ab1dc0SHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(anz,Jptr,lnk);CHKERRQ(ierr);
1864d6ab1dc0SHong Zhang     }
1865d6ab1dc0SHong Zhang 
1866d6ab1dc0SHong Zhang     /* add received col data into lnk */
1867d6ab1dc0SHong Zhang     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1868d6ab1dc0SHong Zhang       if (i == *nextrow[k]) { /* i-th row */
1869d6ab1dc0SHong Zhang         nzi  = *(nextci[k]+1) - *nextci[k];
1870d6ab1dc0SHong Zhang         Jptr = buf_rj[k] + *nextci[k];
1871d6ab1dc0SHong Zhang         ierr = PetscLLCondensedAddSorted_Scalable(nzi,Jptr,lnk);CHKERRQ(ierr);
1872d6ab1dc0SHong Zhang         nextrow[k]++; nextci[k]++;
1873d6ab1dc0SHong Zhang       }
1874d6ab1dc0SHong Zhang     }
1875d6ab1dc0SHong Zhang     nnz = lnk[0];
1876d6ab1dc0SHong Zhang 
1877d6ab1dc0SHong Zhang     /* if free space is not available, make more free space */
1878d6ab1dc0SHong Zhang     if (current_space->local_remaining<nnz) {
1879f91af8c7SBarry Smith       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(nnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
1880d6ab1dc0SHong Zhang       nspacedouble++;
1881d6ab1dc0SHong Zhang     }
1882d6ab1dc0SHong Zhang     /* copy data into free space, then initialize lnk */
1883d6ab1dc0SHong Zhang     ierr = PetscLLCondensedClean_Scalable(nnz,current_space->array,lnk);CHKERRQ(ierr);
1884d6ab1dc0SHong Zhang     ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr);
18852205254eSKarl Rupp 
1886d6ab1dc0SHong Zhang     current_space->array           += nnz;
1887d6ab1dc0SHong Zhang     current_space->local_used      += nnz;
1888d6ab1dc0SHong Zhang     current_space->local_remaining -= nnz;
18892205254eSKarl Rupp 
1890d6ab1dc0SHong Zhang     bi[i+1] = bi[i] + nnz;
1891d6ab1dc0SHong Zhang     if (nnz > rmax) rmax = nnz;
1892d6ab1dc0SHong Zhang   }
1893f671be37SHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1894d6ab1dc0SHong Zhang 
1895854ce69bSBarry Smith   ierr      = PetscMalloc1(bi[pn]+1,&bj);CHKERRQ(ierr);
1896d6ab1dc0SHong Zhang   ierr      = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
1897118727c9SMark F. Adams   afill_tmp = (PetscReal)bi[pn]/(pdti[pn] + poti[pon] + ai[am]+1);
1898d6ab1dc0SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
1899d6ab1dc0SHong Zhang   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr);
19000acc65b4SHong Zhang   ierr = PetscTableDestroy(&ta);CHKERRQ(ierr);
19010acc65b4SHong Zhang 
1902d6ab1dc0SHong Zhang   ierr = MatDestroy(&POt);CHKERRQ(ierr);
1903d6ab1dc0SHong Zhang   ierr = MatDestroy(&PDt);CHKERRQ(ierr);
1904d6ab1dc0SHong Zhang 
1905d6ab1dc0SHong Zhang   /* create symbolic parallel matrix Cmpi - why cannot be assembled in Numeric part   */
1906d6ab1dc0SHong Zhang   /*----------------------------------------------------------------------------------*/
19071795a4d1SJed Brown   ierr = PetscCalloc1(rmax+1,&vals);CHKERRQ(ierr);
1908d6ab1dc0SHong Zhang 
1909d6ab1dc0SHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
1910d6ab1dc0SHong Zhang   ierr = MatSetSizes(Cmpi,pn,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
191133d57670SJed Brown   ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(A->cmap->bs));CHKERRQ(ierr);
1912d6ab1dc0SHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
1913d6ab1dc0SHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
1914d6ab1dc0SHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
1915d6ab1dc0SHong Zhang   ierr = MatSetBlockSize(Cmpi,1);CHKERRQ(ierr);
1916d6ab1dc0SHong Zhang   for (i=0; i<pn; i++) {
1917d6ab1dc0SHong Zhang     row  = i + rstart;
1918d6ab1dc0SHong Zhang     nnz  = bi[i+1] - bi[i];
1919d6ab1dc0SHong Zhang     Jptr = bj + bi[i];
1920d6ab1dc0SHong Zhang     ierr = MatSetValues(Cmpi,1,&row,nnz,Jptr,vals,INSERT_VALUES);CHKERRQ(ierr);
1921d6ab1dc0SHong Zhang   }
1922d6ab1dc0SHong Zhang   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1923d6ab1dc0SHong Zhang   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1924d6ab1dc0SHong Zhang   ierr = PetscFree(vals);CHKERRQ(ierr);
1925d6ab1dc0SHong Zhang 
1926d6ab1dc0SHong Zhang   merge->bi        = bi;
1927d6ab1dc0SHong Zhang   merge->bj        = bj;
1928d6ab1dc0SHong Zhang   merge->coi       = coi;
1929d6ab1dc0SHong Zhang   merge->coj       = coj;
1930d6ab1dc0SHong Zhang   merge->buf_ri    = buf_ri;
1931d6ab1dc0SHong Zhang   merge->buf_rj    = buf_rj;
1932d6ab1dc0SHong Zhang   merge->owners_co = owners_co;
1933d6ab1dc0SHong Zhang 
1934d6ab1dc0SHong Zhang   /* attach the supporting struct to Cmpi for reuse */
1935d6ab1dc0SHong Zhang   c = (Mat_MPIAIJ*)Cmpi->data;
19362205254eSKarl Rupp 
1937d6ab1dc0SHong Zhang   c->ptap     = ptap;
19380298fd71SBarry Smith   ptap->api   = NULL;
19390298fd71SBarry Smith   ptap->apj   = NULL;
1940d6ab1dc0SHong Zhang   ptap->merge = merge;
19410298fd71SBarry Smith   ptap->apa   = NULL;
1942a560ef98SHong Zhang   ptap->destroy   = Cmpi->ops->destroy;
1943a560ef98SHong Zhang   ptap->duplicate = Cmpi->ops->duplicate;
1944a560ef98SHong Zhang 
1945a560ef98SHong Zhang   Cmpi->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ;
1946a560ef98SHong Zhang   Cmpi->ops->destroy                 = MatDestroy_MPIAIJ_PtAP;
1947a560ef98SHong Zhang   Cmpi->ops->duplicate               = MatDuplicate_MPIAIJ_MatPtAP;
1948f96d37fcSHong Zhang 
1949f96d37fcSHong Zhang   *C = Cmpi;
1950f96d37fcSHong Zhang #if defined(PETSC_USE_INFO)
1951f96d37fcSHong Zhang   if (bi[pn] != 0) {
195257622a8eSBarry Smith     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
195357622a8eSBarry Smith     ierr = PetscInfo1(Cmpi,"Use MatTransposeMatMult(A,B,MatReuse,%g,&C) for best performance.\n",(double)afill);CHKERRQ(ierr);
1954f96d37fcSHong Zhang   } else {
1955f96d37fcSHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
1956f96d37fcSHong Zhang   }
1957f96d37fcSHong Zhang #endif
1958f96d37fcSHong Zhang   PetscFunctionReturn(0);
1959f96d37fcSHong Zhang }
1960