xref: /petsc/src/mat/impls/aij/mpi/mpimatmatmult.c (revision c12557a13f0fad7313f160073e352ba05c6ae3a8)
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"
152499ec78SHong Zhang 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 */
202499ec78SHong Zhang 
212499ec78SHong Zhang   PetscFunctionBegin;
222499ec78SHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
230fc8cf34SHong Zhang     ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr);
240fc8cf34SHong Zhang     ierr = PetscOptionsEList("-matmatmult_via","Algorithmic approach","MatMatMult",algTypes,2,algTypes[0],&alg,NULL);CHKERRQ(ierr);
250fc8cf34SHong Zhang     ierr = PetscOptionsEnd();CHKERRQ(ierr);
260fc8cf34SHong Zhang 
273ff4c91cSHong Zhang     ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
280fc8cf34SHong Zhang     switch (alg) {
290fc8cf34SHong Zhang     case 1:
300fc8cf34SHong Zhang       ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(A,B,fill,C);CHKERRQ(ierr);
310fc8cf34SHong Zhang       break;
320fc8cf34SHong Zhang     default:
33b2405163SHong Zhang       ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ(A,B,fill,C);CHKERRQ(ierr);
340fc8cf34SHong Zhang       break;
350fc8cf34SHong Zhang     }
363ff4c91cSHong Zhang     ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
372499ec78SHong Zhang   }
383ff4c91cSHong Zhang   ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
39598bc09dSHong Zhang   ierr = (*(*C)->ops->matmultnumeric)(A,B,*C);CHKERRQ(ierr);
403ff4c91cSHong Zhang   ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
412499ec78SHong Zhang   PetscFunctionReturn(0);
422499ec78SHong Zhang }
432499ec78SHong Zhang 
44f33d1a9aSHong Zhang #undef __FUNCT__
45a1a4e84aSHong Zhang #define __FUNCT__ "MatDestroy_MPIAIJ_MatMatMult"
46a1a4e84aSHong Zhang PetscErrorCode MatDestroy_MPIAIJ_MatMatMult(Mat A)
47598bc09dSHong Zhang {
48598bc09dSHong Zhang   PetscErrorCode ierr;
49598bc09dSHong Zhang   Mat_MPIAIJ     *a    = (Mat_MPIAIJ*)A->data;
50598bc09dSHong Zhang   Mat_PtAPMPI    *ptap = a->ptap;
51598bc09dSHong Zhang 
52598bc09dSHong Zhang   PetscFunctionBegin;
53b7f45c76SHong Zhang   ierr = PetscFree2(ptap->startsj_s,ptap->startsj_r);CHKERRQ(ierr);
54598bc09dSHong Zhang   ierr = PetscFree(ptap->bufa);CHKERRQ(ierr);
55a1a4e84aSHong Zhang   ierr = MatDestroy(&ptap->P_loc);CHKERRQ(ierr);
56a1a4e84aSHong Zhang   ierr = MatDestroy(&ptap->P_oth);CHKERRQ(ierr);
57ab1b013aSHong Zhang   ierr = MatDestroy(&ptap->Pt);CHKERRQ(ierr);
58a1a4e84aSHong Zhang   ierr = PetscFree(ptap->api);CHKERRQ(ierr);
59a1a4e84aSHong Zhang   ierr = PetscFree(ptap->apj);CHKERRQ(ierr);
60598bc09dSHong Zhang   ierr = PetscFree(ptap->apa);CHKERRQ(ierr);
61598bc09dSHong Zhang   ierr = ptap->destroy(A);CHKERRQ(ierr);
62598bc09dSHong Zhang   ierr = PetscFree(ptap);CHKERRQ(ierr);
63598bc09dSHong Zhang   PetscFunctionReturn(0);
64598bc09dSHong Zhang }
65598bc09dSHong Zhang 
662499ec78SHong Zhang #undef __FUNCT__
67a1a4e84aSHong Zhang #define __FUNCT__ "MatDuplicate_MPIAIJ_MatMatMult"
68a1a4e84aSHong Zhang PetscErrorCode MatDuplicate_MPIAIJ_MatMatMult(Mat A, MatDuplicateOption op, Mat *M)
694ae0847bSHong Zhang {
704ae0847bSHong Zhang   PetscErrorCode ierr;
714ae0847bSHong Zhang   Mat_MPIAIJ     *a    = (Mat_MPIAIJ*)A->data;
724ae0847bSHong Zhang   Mat_PtAPMPI    *ptap = a->ptap;
734ae0847bSHong Zhang 
744ae0847bSHong Zhang   PetscFunctionBegin;
754ae0847bSHong Zhang   ierr = (*ptap->duplicate)(A,op,M);CHKERRQ(ierr);
7626fbe8dcSKarl Rupp 
774ae0847bSHong Zhang   (*M)->ops->destroy   = ptap->destroy;   /* = MatDestroy_MPIAIJ, *M doesn't duplicate A's special structure! */
784ae0847bSHong Zhang   (*M)->ops->duplicate = ptap->duplicate; /* = MatDuplicate_MPIAIJ */
794ae0847bSHong Zhang   PetscFunctionReturn(0);
804ae0847bSHong Zhang }
814ae0847bSHong Zhang 
82*c12557a1SHong Zhang /* compute apa = A[i,:]*P */
83*c12557a1SHong Zhang #define AProw(i,ad,ao,p_loc,p_oth,apa) \
84*c12557a1SHong Zhang {\
85*c12557a1SHong Zhang   PetscInt    _anz,_pnz,_j,_k,*_ai,*_aj,_row,*_pi,*_pj;      \
86*c12557a1SHong Zhang   PetscScalar *_aa,_valtmp,*_pa;                             \
87*c12557a1SHong Zhang   /* diagonal portion of A */\
88*c12557a1SHong Zhang   _ai  = ad->i;\
89*c12557a1SHong Zhang   _anz = _ai[i+1] - _ai[i];\
90*c12557a1SHong Zhang   _aj  = ad->j + _ai[i];\
91*c12557a1SHong Zhang   _aa  = ad->a + _ai[i];\
92*c12557a1SHong Zhang   for (_j=0; _j<_anz; _j++) {\
93*c12557a1SHong Zhang     _row = _aj[_j]; \
94*c12557a1SHong Zhang     _pi  = p_loc->i;                                 \
95*c12557a1SHong Zhang     _pnz = _pi[_row+1] - _pi[_row];         \
96*c12557a1SHong Zhang     _pj  = p_loc->j + _pi[_row];                 \
97*c12557a1SHong Zhang     _pa  = p_loc->a + _pi[_row];                 \
98*c12557a1SHong Zhang     /* perform dense axpy */                    \
99*c12557a1SHong Zhang     valtmp = _aa[_j];                           \
100*c12557a1SHong Zhang     for (_k=0; _k<_pnz; _k++) {                    \
101*c12557a1SHong Zhang       apa[_pj[_k]] += valtmp*_pa[_k];               \
102*c12557a1SHong Zhang     }                                           \
103*c12557a1SHong Zhang     PetscLogFlops(2.0*_pnz);                    \
104*c12557a1SHong Zhang   }                                             \
105*c12557a1SHong Zhang   /* off-diagonal portion of A */               \
106*c12557a1SHong Zhang   _ai  = ao->i;\
107*c12557a1SHong Zhang   _anz = _ai[i+1] - _ai[i];                     \
108*c12557a1SHong Zhang   _aj  = ao->j + _ai[i];                         \
109*c12557a1SHong Zhang   _aa  = ao->a + _ai[i];                         \
110*c12557a1SHong Zhang   for (_j=0; _j<_anz; _j++) {                      \
111*c12557a1SHong Zhang     _row = _aj[_j];    \
112*c12557a1SHong Zhang     _pi  = p_oth->i;                         \
113*c12557a1SHong Zhang     _pnz = _pi[_row+1] - _pi[_row];          \
114*c12557a1SHong Zhang     _pj  = p_oth->j + _pi[_row];                  \
115*c12557a1SHong Zhang     _pa  = p_oth->a + _pi[_row];                  \
116*c12557a1SHong Zhang     /* perform dense axpy */                     \
117*c12557a1SHong Zhang     _valtmp = _aa[_j];                             \
118*c12557a1SHong Zhang     for (_k=0; _k<_pnz; _k++) {                     \
119*c12557a1SHong Zhang       apa[_pj[_k]] += _valtmp*_pa[_k];                \
120*c12557a1SHong Zhang     }                                            \
121*c12557a1SHong Zhang     PetscLogFlops(2.0*_pnz);                     \
122*c12557a1SHong Zhang   }                                              \
123*c12557a1SHong Zhang }
124*c12557a1SHong Zhang 
1254ae0847bSHong Zhang #undef __FUNCT__
126f8efd71cSHong Zhang #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable"
127f8efd71cSHong Zhang PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable(Mat A,Mat P,Mat C)
128598bc09dSHong Zhang {
129598bc09dSHong Zhang   PetscErrorCode ierr;
1304ae0847bSHong Zhang   Mat_MPIAIJ     *a  =(Mat_MPIAIJ*)A->data,*c=(Mat_MPIAIJ*)C->data;
131598bc09dSHong Zhang   Mat_SeqAIJ     *ad =(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
132598bc09dSHong Zhang   Mat_SeqAIJ     *cd =(Mat_SeqAIJ*)(c->A)->data,*co=(Mat_SeqAIJ*)(c->B)->data;
133*c12557a1SHong Zhang   PetscScalar    *cda=cd->a,*coa=co->a;
134598bc09dSHong Zhang   Mat_SeqAIJ     *p_loc,*p_oth;
135*c12557a1SHong Zhang   PetscScalar    *apa,valtmp,*ca;
136*c12557a1SHong Zhang   PetscInt       cm   =C->rmap->n;
137598bc09dSHong Zhang   Mat_PtAPMPI    *ptap=c->ptap;
138*c12557a1SHong Zhang   PetscInt       *api,*apj,*apJ,i,k;
13929825010SHong Zhang   PetscInt       cstart=C->cmap->rstart;
140598bc09dSHong Zhang   PetscInt       cdnz,conz,k0,k1;
1419467f45fSHong Zhang   MPI_Comm       comm;
1429467f45fSHong Zhang   PetscMPIInt    size;
143598bc09dSHong Zhang 
144598bc09dSHong Zhang   PetscFunctionBegin;
1459467f45fSHong Zhang   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
1469467f45fSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1479467f45fSHong Zhang 
148a1a4e84aSHong Zhang   /* 1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
149598bc09dSHong Zhang   /*-----------------------------------------------------*/
150a1a4e84aSHong Zhang   /* update numerical values of P_oth and P_loc */
151b7f45c76SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
152a1a4e84aSHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
153598bc09dSHong Zhang 
154598bc09dSHong Zhang   /* 2) compute numeric C_loc = A_loc*P = Ad*P_loc + Ao*P_oth */
155598bc09dSHong Zhang   /*----------------------------------------------------------*/
156598bc09dSHong Zhang   /* get data from symbolic products */
157a1a4e84aSHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
158*c12557a1SHong Zhang   p_oth = NULL;
1599467f45fSHong Zhang   if (size >1) {
1609467f45fSHong Zhang     p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
1619467f45fSHong Zhang   }
162598bc09dSHong Zhang 
163598bc09dSHong Zhang   /* get apa for storing dense row A[i,:]*P */
164598bc09dSHong Zhang   apa = ptap->apa;
165598bc09dSHong Zhang 
16629825010SHong Zhang   api = ptap->api;
16729825010SHong Zhang   apj = ptap->apj;
168598bc09dSHong Zhang   for (i=0; i<cm; i++) {
169*c12557a1SHong Zhang     /* compute apa = A[i,:]*P */
170*c12557a1SHong Zhang     AProw(i,ad,ao,p_loc,p_oth,apa);
171598bc09dSHong Zhang 
172598bc09dSHong Zhang     /* set values in C */
173598bc09dSHong Zhang     apJ  = apj + api[i];
174598bc09dSHong Zhang     cdnz = cd->i[i+1] - cd->i[i];
175598bc09dSHong Zhang     conz = co->i[i+1] - co->i[i];
176598bc09dSHong Zhang 
177598bc09dSHong Zhang     /* 1st off-diagoanl part of C */
178598bc09dSHong Zhang     ca = coa + co->i[i];
179598bc09dSHong Zhang     k  = 0;
180598bc09dSHong Zhang     for (k0=0; k0<conz; k0++) {
181598bc09dSHong Zhang       if (apJ[k] >= cstart) break;
182598bc09dSHong Zhang       ca[k0]      = apa[apJ[k]];
183598bc09dSHong Zhang       apa[apJ[k]] = 0.0;
184598bc09dSHong Zhang       k++;
185598bc09dSHong Zhang     }
186598bc09dSHong Zhang 
187598bc09dSHong Zhang     /* diagonal part of C */
188598bc09dSHong Zhang     ca = cda + cd->i[i];
189598bc09dSHong Zhang     for (k1=0; k1<cdnz; k1++) {
190598bc09dSHong Zhang       ca[k1]      = apa[apJ[k]];
191598bc09dSHong Zhang       apa[apJ[k]] = 0.0;
192598bc09dSHong Zhang       k++;
193598bc09dSHong Zhang     }
194598bc09dSHong Zhang 
195598bc09dSHong Zhang     /* 2nd off-diagoanl part of C */
196598bc09dSHong Zhang     ca = coa + co->i[i];
197598bc09dSHong Zhang     for (; k0<conz; k0++) {
198598bc09dSHong Zhang       ca[k0]      = apa[apJ[k]];
199598bc09dSHong Zhang       apa[apJ[k]] = 0.0;
200598bc09dSHong Zhang       k++;
201598bc09dSHong Zhang     }
202598bc09dSHong Zhang   }
203598bc09dSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
204598bc09dSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
205598bc09dSHong Zhang   PetscFunctionReturn(0);
206598bc09dSHong Zhang }
207598bc09dSHong Zhang 
208598bc09dSHong Zhang #undef __FUNCT__
2090fc8cf34SHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable"
2100fc8cf34SHong Zhang PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat A,Mat P,PetscReal fill,Mat *C)
211598bc09dSHong Zhang {
212598bc09dSHong Zhang   PetscErrorCode     ierr;
213ce94432eSBarry Smith   MPI_Comm           comm;
2149467f45fSHong Zhang   PetscMPIInt        size;
215bfcd1a73SHong Zhang   Mat                Cmpi;
216598bc09dSHong Zhang   Mat_PtAPMPI        *ptap;
2170298fd71SBarry Smith   PetscFreeSpaceList free_space=NULL,current_space=NULL;
2184ae0847bSHong Zhang   Mat_MPIAIJ         *a        =(Mat_MPIAIJ*)A->data,*c;
219bfcd1a73SHong Zhang   Mat_SeqAIJ         *ad       =(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_loc,*p_oth;
2204ae0847bSHong Zhang   PetscInt           *pi_loc,*pj_loc,*pi_oth,*pj_oth,*dnz,*onz;
2214ae0847bSHong Zhang   PetscInt           *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,rstart=A->rmap->rstart;
222d14fa243SHong Zhang   PetscInt           *lnk,i,pnz,row,*api,*apj,*Jptr,apnz,nspacedouble=0,j,nzi;
223bfcd1a73SHong Zhang   PetscInt           am=A->rmap->n,pN=P->cmap->N,pn=P->cmap->n,pm=P->rmap->n;
224598bc09dSHong Zhang   PetscBT            lnkbt;
225598bc09dSHong Zhang   PetscScalar        *apa;
226bfcd1a73SHong Zhang   PetscReal          afill;
227f84c536bSHong Zhang   PetscInt           nlnk_max,armax,prmax;
228598bc09dSHong Zhang 
229598bc09dSHong Zhang   PetscFunctionBegin;
230ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
2319467f45fSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
2329467f45fSHong Zhang 
233a1a4e84aSHong Zhang   if (A->cmap->rstart != P->rmap->rstart || A->cmap->rend != P->rmap->rend) {
234cf1a0672SHong Zhang     SETERRQ4(comm,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, (%D, %D) != (%D,%D)",A->cmap->rstart,A->cmap->rend,P->rmap->rstart,P->rmap->rend);
235a1a4e84aSHong Zhang   }
236152983bfSHong Zhang 
237598bc09dSHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
238b00a9115SJed Brown   ierr = PetscNew(&ptap);CHKERRQ(ierr);
239598bc09dSHong Zhang 
240598bc09dSHong Zhang   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
241b7f45c76SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
24219f0e57aSHong Zhang 
243598bc09dSHong Zhang   /* get P_loc by taking all local rows of P */
244a1a4e84aSHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
245598bc09dSHong Zhang 
246a1a4e84aSHong Zhang   p_loc  = (Mat_SeqAIJ*)(ptap->P_loc)->data;
247598bc09dSHong Zhang   pi_loc = p_loc->i; pj_loc = p_loc->j;
2489467f45fSHong Zhang   if (size > 1) {
2499467f45fSHong Zhang     p_oth  = (Mat_SeqAIJ*)(ptap->P_oth)->data;
250598bc09dSHong Zhang     pi_oth = p_oth->i; pj_oth = p_oth->j;
25120e3dc75SHong Zhang   } else {
25220e3dc75SHong Zhang     p_oth = NULL;
25320e3dc75SHong Zhang     pi_oth = NULL; pj_oth = NULL;
2549467f45fSHong Zhang   }
255598bc09dSHong Zhang 
256598bc09dSHong Zhang   /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */
257598bc09dSHong Zhang   /*-------------------------------------------------------------------*/
258854ce69bSBarry Smith   ierr      = PetscMalloc1(am+2,&api);CHKERRQ(ierr);
259a1a4e84aSHong Zhang   ptap->api = api;
260598bc09dSHong Zhang   api[0]    = 0;
261598bc09dSHong Zhang 
262598bc09dSHong Zhang   /* create and initialize a linked list */
263f84c536bSHong Zhang   armax    = ad->rmax+ao->rmax;
2649467f45fSHong Zhang   if (size >1) {
265f84c536bSHong Zhang     prmax    = PetscMax(p_loc->rmax,p_oth->rmax);
2669467f45fSHong Zhang   } else {
2679467f45fSHong Zhang     prmax = p_loc->rmax;
2689467f45fSHong Zhang   }
269f84c536bSHong Zhang   nlnk_max = armax*prmax;
270f84c536bSHong Zhang   if (!nlnk_max || nlnk_max > pN) nlnk_max = pN;
2710d3134e9SHong Zhang   ierr = PetscLLCondensedCreate(nlnk_max,pN,&lnk,&lnkbt);CHKERRQ(ierr);
272598bc09dSHong Zhang 
273bfcd1a73SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */
274bfcd1a73SHong Zhang   ierr = PetscFreeSpaceGet((PetscInt)(fill*(adi[am]+aoi[am]+pi_loc[pm])),&free_space);CHKERRQ(ierr);
2752205254eSKarl Rupp 
276598bc09dSHong Zhang   current_space = free_space;
277598bc09dSHong Zhang 
278598bc09dSHong Zhang   ierr = MatPreallocateInitialize(comm,am,pn,dnz,onz);CHKERRQ(ierr);
279598bc09dSHong Zhang   for (i=0; i<am; i++) {
280598bc09dSHong Zhang     /* diagonal portion of A */
281598bc09dSHong Zhang     nzi = adi[i+1] - adi[i];
282598bc09dSHong Zhang     for (j=0; j<nzi; j++) {
283598bc09dSHong Zhang       row  = *adj++;
284598bc09dSHong Zhang       pnz  = pi_loc[row+1] - pi_loc[row];
285598bc09dSHong Zhang       Jptr = pj_loc + pi_loc[row];
286598bc09dSHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
2871fbbb359SHong Zhang       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
288598bc09dSHong Zhang     }
289598bc09dSHong Zhang     /* off-diagonal portion of A */
290598bc09dSHong Zhang     nzi = aoi[i+1] - aoi[i];
291598bc09dSHong Zhang     for (j=0; j<nzi; j++) {
292598bc09dSHong Zhang       row  = *aoj++;
293598bc09dSHong Zhang       pnz  = pi_oth[row+1] - pi_oth[row];
294598bc09dSHong Zhang       Jptr = pj_oth + pi_oth[row];
2951fbbb359SHong Zhang       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
296598bc09dSHong Zhang     }
297598bc09dSHong Zhang 
298d14fa243SHong Zhang     apnz     = lnk[0];
299598bc09dSHong Zhang     api[i+1] = api[i] + apnz;
300598bc09dSHong Zhang 
301598bc09dSHong Zhang     /* if free space is not available, double the total space in the list */
302598bc09dSHong Zhang     if (current_space->local_remaining<apnz) {
303598bc09dSHong Zhang       ierr = PetscFreeSpaceGet(apnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
304598bc09dSHong Zhang       nspacedouble++;
305598bc09dSHong Zhang     }
306598bc09dSHong Zhang 
307598bc09dSHong Zhang     /* Copy data into free space, then initialize lnk */
3081fbbb359SHong Zhang     ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
309598bc09dSHong Zhang     ierr = MatPreallocateSet(i+rstart,apnz,current_space->array,dnz,onz);CHKERRQ(ierr);
3102205254eSKarl Rupp 
311598bc09dSHong Zhang     current_space->array           += apnz;
312598bc09dSHong Zhang     current_space->local_used      += apnz;
313598bc09dSHong Zhang     current_space->local_remaining -= apnz;
314598bc09dSHong Zhang   }
315598bc09dSHong Zhang 
316598bc09dSHong Zhang   /* Allocate space for apj, initialize apj, and */
317598bc09dSHong Zhang   /* destroy list of free space and other temporary array(s) */
318854ce69bSBarry Smith   ierr = PetscMalloc1(api[am]+1,&ptap->apj);CHKERRQ(ierr);
319a1a4e84aSHong Zhang   apj  = ptap->apj;
320a1a4e84aSHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,ptap->apj);CHKERRQ(ierr);
321598bc09dSHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
322598bc09dSHong Zhang 
323f84c536bSHong Zhang   /* malloc apa to store dense row A[i,:]*P */
3241795a4d1SJed Brown   ierr = PetscCalloc1(pN,&apa);CHKERRQ(ierr);
3252205254eSKarl Rupp 
326f84c536bSHong Zhang   ptap->apa = apa;
327f84c536bSHong Zhang 
328bfcd1a73SHong Zhang   /* create and assemble symbolic parallel matrix Cmpi */
329598bc09dSHong Zhang   /*----------------------------------------------------*/
330bfcd1a73SHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
331bfcd1a73SHong Zhang   ierr = MatSetSizes(Cmpi,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
33233d57670SJed Brown   ierr = MatSetBlockSizesFromMats(Cmpi,A,P);CHKERRQ(ierr);
333a2f3521dSMark F. Adams 
334bfcd1a73SHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
335bfcd1a73SHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
336598bc09dSHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
337598bc09dSHong Zhang   for (i=0; i<am; i++) {
338598bc09dSHong Zhang     row  = i + rstart;
339598bc09dSHong Zhang     apnz = api[i+1] - api[i];
340bfcd1a73SHong Zhang     ierr = MatSetValues(Cmpi,1,&row,apnz,apj,apa,INSERT_VALUES);CHKERRQ(ierr);
341598bc09dSHong Zhang     apj += apnz;
342598bc09dSHong Zhang   }
343bfcd1a73SHong Zhang   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
344bfcd1a73SHong Zhang   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
345598bc09dSHong Zhang 
346bfcd1a73SHong Zhang   ptap->destroy        = Cmpi->ops->destroy;
347bfcd1a73SHong Zhang   ptap->duplicate      = Cmpi->ops->duplicate;
348f8efd71cSHong Zhang   Cmpi->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable;
349bfcd1a73SHong Zhang   Cmpi->ops->destroy   = MatDestroy_MPIAIJ_MatMatMult;
350bfcd1a73SHong Zhang   Cmpi->ops->duplicate = MatDuplicate_MPIAIJ_MatMatMult;
351598bc09dSHong Zhang 
352bfcd1a73SHong Zhang   /* attach the supporting struct to Cmpi for reuse */
353bfcd1a73SHong Zhang   c       = (Mat_MPIAIJ*)Cmpi->data;
354598bc09dSHong Zhang   c->ptap = ptap;
355598bc09dSHong Zhang 
356bfcd1a73SHong Zhang   *C = Cmpi;
357bfcd1a73SHong Zhang 
358bfcd1a73SHong Zhang   /* set MatInfo */
359118727c9SMark F. Adams   afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1) + 1.e-5;
360bfcd1a73SHong Zhang   if (afill < 1.0) afill = 1.0;
361bfcd1a73SHong Zhang   Cmpi->info.mallocs           = nspacedouble;
362bfcd1a73SHong Zhang   Cmpi->info.fill_ratio_given  = fill;
363bfcd1a73SHong Zhang   Cmpi->info.fill_ratio_needed = afill;
364bfcd1a73SHong Zhang 
365bfcd1a73SHong Zhang #if defined(PETSC_USE_INFO)
366bfcd1a73SHong Zhang   if (api[am]) {
36757622a8eSBarry Smith     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
36857622a8eSBarry Smith     ierr = PetscInfo1(Cmpi,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr);
369bfcd1a73SHong Zhang   } else {
370bfcd1a73SHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
371bfcd1a73SHong Zhang   }
372bfcd1a73SHong Zhang #endif
373598bc09dSHong Zhang   PetscFunctionReturn(0);
374598bc09dSHong Zhang }
375598bc09dSHong Zhang 
3769123193aSHong Zhang #undef __FUNCT__
3779123193aSHong Zhang #define __FUNCT__ "MatMatMult_MPIAIJ_MPIDense"
3789123193aSHong Zhang PetscErrorCode MatMatMult_MPIAIJ_MPIDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
3799123193aSHong Zhang {
3809123193aSHong Zhang   PetscErrorCode ierr;
3819123193aSHong Zhang 
3829123193aSHong Zhang   PetscFunctionBegin;
3839123193aSHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
3843ff4c91cSHong Zhang     ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
3859123193aSHong Zhang     ierr = MatMatMultSymbolic_MPIAIJ_MPIDense(A,B,fill,C);CHKERRQ(ierr);
3863ff4c91cSHong Zhang     ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
3879123193aSHong Zhang   }
3883ff4c91cSHong Zhang   ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
3899123193aSHong Zhang   ierr = MatMatMultNumeric_MPIAIJ_MPIDense(A,B,*C);CHKERRQ(ierr);
3903ff4c91cSHong Zhang   ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
3919123193aSHong Zhang   PetscFunctionReturn(0);
3929123193aSHong Zhang }
3939123193aSHong Zhang 
3944324174eSBarry Smith typedef struct {
3954324174eSBarry Smith   Mat         workB;
3964324174eSBarry Smith   PetscScalar *rvalues,*svalues;
3974324174eSBarry Smith   MPI_Request *rwaits,*swaits;
3984324174eSBarry Smith } MPIAIJ_MPIDense;
3994324174eSBarry Smith 
4007af9e4fdSHong Zhang #undef __FUNCT__
40196b95a6bSBarry Smith #define __FUNCT__ "MatMPIAIJ_MPIDenseDestroy"
40296b95a6bSBarry Smith PetscErrorCode MatMPIAIJ_MPIDenseDestroy(void *ctx)
4034324174eSBarry Smith {
4044324174eSBarry Smith   MPIAIJ_MPIDense *contents = (MPIAIJ_MPIDense*) ctx;
4054324174eSBarry Smith   PetscErrorCode  ierr;
4064324174eSBarry Smith 
4074324174eSBarry Smith   PetscFunctionBegin;
4086bf464f9SBarry Smith   ierr = MatDestroy(&contents->workB);CHKERRQ(ierr);
4094324174eSBarry Smith   ierr = PetscFree4(contents->rvalues,contents->svalues,contents->rwaits,contents->swaits);CHKERRQ(ierr);
4104324174eSBarry Smith   ierr = PetscFree(contents);CHKERRQ(ierr);
4114324174eSBarry Smith   PetscFunctionReturn(0);
4124324174eSBarry Smith }
4134324174eSBarry Smith 
4149123193aSHong Zhang #undef __FUNCT__
415fd4e9aacSBarry Smith #define __FUNCT__ "MatMatMultNumeric_MPIDense"
416fd4e9aacSBarry Smith /*
417fd4e9aacSBarry Smith     This is a "dummy function" that handles the case where matrix C was created as a dense matrix
418fd4e9aacSBarry Smith   directly by the user and passed to MatMatMult() with the MAT_REUSE_MATRIX option
419fd4e9aacSBarry Smith 
420fd4e9aacSBarry Smith   It is the same as MatMatMultSymbolic_MPIAIJ_MPIDense() except does not create C
421fd4e9aacSBarry Smith */
422fd4e9aacSBarry Smith PetscErrorCode MatMatMultNumeric_MPIDense(Mat A,Mat B,Mat C)
423fd4e9aacSBarry Smith {
424fd4e9aacSBarry Smith   PetscErrorCode         ierr;
425fd4e9aacSBarry Smith   PetscBool              flg;
426fd4e9aacSBarry Smith   Mat_MPIAIJ             *aij = (Mat_MPIAIJ*) A->data;
427fd4e9aacSBarry Smith   PetscInt               nz   = aij->B->cmap->n;
428fd4e9aacSBarry Smith   PetscContainer         container;
429fd4e9aacSBarry Smith   MPIAIJ_MPIDense        *contents;
430fd4e9aacSBarry Smith   VecScatter             ctx   = aij->Mvctx;
431fd4e9aacSBarry Smith   VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata;
432fd4e9aacSBarry Smith   VecScatter_MPI_General *to   = (VecScatter_MPI_General*) ctx->todata;
433fd4e9aacSBarry Smith 
434fd4e9aacSBarry Smith   PetscFunctionBegin;
435fd4e9aacSBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)B,MATMPIDENSE,&flg);CHKERRQ(ierr);
436fd4e9aacSBarry Smith   if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Second matrix must be mpidense");
437fd4e9aacSBarry Smith 
438fd4e9aacSBarry Smith   /* Handle case where where user provided the final C matrix rather than calling MatMatMult() with MAT_INITIAL_MATRIX*/
439fd4e9aacSBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&flg);CHKERRQ(ierr);
440fd4e9aacSBarry Smith   if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"First matrix must be MPIAIJ");
441fd4e9aacSBarry Smith 
442fd4e9aacSBarry Smith   C->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIDense;
443fd4e9aacSBarry Smith 
444fd4e9aacSBarry Smith   ierr = PetscNew(&contents);CHKERRQ(ierr);
445fd4e9aacSBarry Smith   /* Create work matrix used to store off processor rows of B needed for local product */
446fd4e9aacSBarry Smith   ierr = MatCreateSeqDense(PETSC_COMM_SELF,nz,B->cmap->N,NULL,&contents->workB);CHKERRQ(ierr);
447fd4e9aacSBarry Smith   /* Create work arrays needed */
448fd4e9aacSBarry Smith   ierr = PetscMalloc4(B->cmap->N*from->starts[from->n],&contents->rvalues,
449fd4e9aacSBarry Smith                       B->cmap->N*to->starts[to->n],&contents->svalues,
450fd4e9aacSBarry Smith                       from->n,&contents->rwaits,
451fd4e9aacSBarry Smith                       to->n,&contents->swaits);CHKERRQ(ierr);
452fd4e9aacSBarry Smith 
453fd4e9aacSBarry Smith   ierr = PetscContainerCreate(PetscObjectComm((PetscObject)A),&container);CHKERRQ(ierr);
454fd4e9aacSBarry Smith   ierr = PetscContainerSetPointer(container,contents);CHKERRQ(ierr);
455fd4e9aacSBarry Smith   ierr = PetscContainerSetUserDestroy(container,MatMPIAIJ_MPIDenseDestroy);CHKERRQ(ierr);
456fd4e9aacSBarry Smith   ierr = PetscObjectCompose((PetscObject)C,"workB",(PetscObject)container);CHKERRQ(ierr);
457fd4e9aacSBarry Smith   ierr = PetscContainerDestroy(&container);CHKERRQ(ierr);
458fd4e9aacSBarry Smith 
459fd4e9aacSBarry Smith   ierr = (*C->ops->matmultnumeric)(A,B,C);CHKERRQ(ierr);
460fd4e9aacSBarry Smith   PetscFunctionReturn(0);
461fd4e9aacSBarry Smith }
462fd4e9aacSBarry Smith 
463fd4e9aacSBarry Smith #undef __FUNCT__
4649123193aSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIDense"
4659123193aSHong Zhang PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIDense(Mat A,Mat B,PetscReal fill,Mat *C)
4669123193aSHong Zhang {
4679123193aSHong Zhang   PetscErrorCode         ierr;
468f32d5d43SBarry Smith   Mat_MPIAIJ             *aij = (Mat_MPIAIJ*) A->data;
469d0f46423SBarry Smith   PetscInt               nz   = aij->B->cmap->n;
470bf0cc555SLisandro Dalcin   PetscContainer         container;
4714324174eSBarry Smith   MPIAIJ_MPIDense        *contents;
4724324174eSBarry Smith   VecScatter             ctx   = aij->Mvctx;
4734324174eSBarry Smith   VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata;
4744324174eSBarry Smith   VecScatter_MPI_General *to   = (VecScatter_MPI_General*) ctx->todata;
475d0f46423SBarry Smith   PetscInt               m     = A->rmap->n,n=B->cmap->n;
4769123193aSHong Zhang 
4779123193aSHong Zhang   PetscFunctionBegin;
478ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)B),C);CHKERRQ(ierr);
479d0f46423SBarry Smith   ierr = MatSetSizes(*C,m,n,A->rmap->N,B->cmap->N);CHKERRQ(ierr);
48033d57670SJed Brown   ierr = MatSetBlockSizesFromMats(*C,A,B);CHKERRQ(ierr);
481cb2480eaSBarry Smith   ierr = MatSetType(*C,MATMPIDENSE);CHKERRQ(ierr);
4820298fd71SBarry Smith   ierr = MatMPIDenseSetPreallocation(*C,NULL);CHKERRQ(ierr);
483cb2480eaSBarry Smith   ierr = MatAssemblyBegin(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
484cb2480eaSBarry Smith   ierr = MatAssemblyEnd(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
4852205254eSKarl Rupp 
486d73949e8SHong Zhang   (*C)->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIDense;
487f32d5d43SBarry Smith 
488b00a9115SJed Brown   ierr = PetscNew(&contents);CHKERRQ(ierr);
489f32d5d43SBarry Smith   /* Create work matrix used to store off processor rows of B needed for local product */
4900298fd71SBarry Smith   ierr = MatCreateSeqDense(PETSC_COMM_SELF,nz,B->cmap->N,NULL,&contents->workB);CHKERRQ(ierr);
4914324174eSBarry Smith   /* Create work arrays needed */
492dcca6d9dSJed Brown   ierr = PetscMalloc4(B->cmap->N*from->starts[from->n],&contents->rvalues,
493dcca6d9dSJed Brown                       B->cmap->N*to->starts[to->n],&contents->svalues,
494dcca6d9dSJed Brown                       from->n,&contents->rwaits,
495dcca6d9dSJed Brown                       to->n,&contents->swaits);CHKERRQ(ierr);
4964324174eSBarry Smith 
497ce94432eSBarry Smith   ierr = PetscContainerCreate(PetscObjectComm((PetscObject)A),&container);CHKERRQ(ierr);
498bf0cc555SLisandro Dalcin   ierr = PetscContainerSetPointer(container,contents);CHKERRQ(ierr);
49996b95a6bSBarry Smith   ierr = PetscContainerSetUserDestroy(container,MatMPIAIJ_MPIDenseDestroy);CHKERRQ(ierr);
500bf0cc555SLisandro Dalcin   ierr = PetscObjectCompose((PetscObject)(*C),"workB",(PetscObject)container);CHKERRQ(ierr);
501bf0cc555SLisandro Dalcin   ierr = PetscContainerDestroy(&container);CHKERRQ(ierr);
5029123193aSHong Zhang   PetscFunctionReturn(0);
5039123193aSHong Zhang }
5049123193aSHong Zhang 
5057af9e4fdSHong Zhang #undef __FUNCT__
5067af9e4fdSHong Zhang #define __FUNCT__ "MatMPIDenseScatter"
507f32d5d43SBarry Smith /*
5082636ff26SBarry Smith     Performs an efficient scatter on the rows of B needed by this process; this is
5092636ff26SBarry Smith     a modification of the VecScatterBegin_() routines.
510f32d5d43SBarry Smith */
5114324174eSBarry Smith PetscErrorCode MatMPIDenseScatter(Mat A,Mat B,Mat C,Mat *outworkB)
512f32d5d43SBarry Smith {
513f32d5d43SBarry Smith   Mat_MPIAIJ             *aij = (Mat_MPIAIJ*)A->data;
514f32d5d43SBarry Smith   PetscErrorCode         ierr;
515f32d5d43SBarry Smith   PetscScalar            *b,*w,*svalues,*rvalues;
516f32d5d43SBarry Smith   VecScatter             ctx   = aij->Mvctx;
517f32d5d43SBarry Smith   VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata;
518f32d5d43SBarry Smith   VecScatter_MPI_General *to   = (VecScatter_MPI_General*) ctx->todata;
519f32d5d43SBarry Smith   PetscInt               i,j,k;
520aa5bb8c0SSatish Balay   PetscInt               *sindices,*sstarts,*rindices,*rstarts;
521aa5bb8c0SSatish Balay   PetscMPIInt            *sprocs,*rprocs,nrecvs;
522f32d5d43SBarry Smith   MPI_Request            *swaits,*rwaits;
523ce94432eSBarry Smith   MPI_Comm               comm;
524d0f46423SBarry Smith   PetscMPIInt            tag  = ((PetscObject)ctx)->tag,ncols = B->cmap->N, nrows = aij->B->cmap->n,imdex,nrowsB = B->rmap->n;
525f32d5d43SBarry Smith   MPI_Status             status;
5264324174eSBarry Smith   MPIAIJ_MPIDense        *contents;
527bf0cc555SLisandro Dalcin   PetscContainer         container;
5284324174eSBarry Smith   Mat                    workB;
529f32d5d43SBarry Smith 
530f32d5d43SBarry Smith   PetscFunctionBegin;
531ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
532bf0cc555SLisandro Dalcin   ierr = PetscObjectQuery((PetscObject)C,"workB",(PetscObject*)&container);CHKERRQ(ierr);
53329825010SHong Zhang   if (!container) SETERRQ(comm,PETSC_ERR_PLIB,"Container does not exist");
534bf0cc555SLisandro Dalcin   ierr = PetscContainerGetPointer(container,(void**)&contents);CHKERRQ(ierr);
5354324174eSBarry Smith 
5364324174eSBarry Smith   workB = *outworkB = contents->workB;
537cf1a0672SHong 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);
538f32d5d43SBarry Smith   sindices = to->indices;
539f32d5d43SBarry Smith   sstarts  = to->starts;
540f32d5d43SBarry Smith   sprocs   = to->procs;
5414324174eSBarry Smith   swaits   = contents->swaits;
5424324174eSBarry Smith   svalues  = contents->svalues;
543f32d5d43SBarry Smith 
544f32d5d43SBarry Smith   rindices = from->indices;
545f32d5d43SBarry Smith   rstarts  = from->starts;
546f32d5d43SBarry Smith   rprocs   = from->procs;
5474324174eSBarry Smith   rwaits   = contents->rwaits;
5484324174eSBarry Smith   rvalues  = contents->rvalues;
549f32d5d43SBarry Smith 
5508c778c55SBarry Smith   ierr = MatDenseGetArray(B,&b);CHKERRQ(ierr);
5518c778c55SBarry Smith   ierr = MatDenseGetArray(workB,&w);CHKERRQ(ierr);
552f32d5d43SBarry Smith 
553f32d5d43SBarry Smith   for (i=0; i<from->n; i++) {
554f32d5d43SBarry Smith     ierr = MPI_Irecv(rvalues+ncols*rstarts[i],ncols*(rstarts[i+1]-rstarts[i]),MPIU_SCALAR,rprocs[i],tag,comm,rwaits+i);CHKERRQ(ierr);
555f32d5d43SBarry Smith   }
556f32d5d43SBarry Smith 
557f32d5d43SBarry Smith   for (i=0; i<to->n; i++) {
558f32d5d43SBarry Smith     /* pack a message at a time */
559f32d5d43SBarry Smith     for (j=0; j<sstarts[i+1]-sstarts[i]; j++) {
560f32d5d43SBarry Smith       for (k=0; k<ncols; k++) {
561f32d5d43SBarry Smith         svalues[ncols*(sstarts[i] + j) + k] = b[sindices[sstarts[i]+j] + nrowsB*k];
5622636ff26SBarry Smith       }
5632636ff26SBarry Smith     }
564f32d5d43SBarry Smith     ierr = MPI_Isend(svalues+ncols*sstarts[i],ncols*(sstarts[i+1]-sstarts[i]),MPIU_SCALAR,sprocs[i],tag,comm,swaits+i);CHKERRQ(ierr);
565f32d5d43SBarry Smith   }
5662636ff26SBarry Smith 
567f32d5d43SBarry Smith   nrecvs = from->n;
568f32d5d43SBarry Smith   while (nrecvs) {
569f32d5d43SBarry Smith     ierr = MPI_Waitany(from->n,rwaits,&imdex,&status);CHKERRQ(ierr);
570f32d5d43SBarry Smith     nrecvs--;
571f32d5d43SBarry Smith     /* unpack a message at a time */
572f32d5d43SBarry Smith     for (j=0; j<rstarts[imdex+1]-rstarts[imdex]; j++) {
573f32d5d43SBarry Smith       for (k=0; k<ncols; k++) {
574f32d5d43SBarry Smith         w[rindices[rstarts[imdex]+j] + nrows*k] = rvalues[ncols*(rstarts[imdex] + j) + k];
5752636ff26SBarry Smith       }
5762636ff26SBarry Smith     }
577f32d5d43SBarry Smith   }
578cb9801acSJed Brown   if (to->n) {ierr = MPI_Waitall(to->n,swaits,to->sstatus);CHKERRQ(ierr);}
579f32d5d43SBarry Smith 
5808c778c55SBarry Smith   ierr = MatDenseRestoreArray(B,&b);CHKERRQ(ierr);
5818c778c55SBarry Smith   ierr = MatDenseRestoreArray(workB,&w);CHKERRQ(ierr);
582f32d5d43SBarry Smith   ierr = MatAssemblyBegin(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
583f32d5d43SBarry Smith   ierr = MatAssemblyEnd(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
584f32d5d43SBarry Smith   PetscFunctionReturn(0);
585f32d5d43SBarry Smith }
586f32d5d43SBarry Smith extern PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat,Mat,Mat);
587f32d5d43SBarry Smith 
5889123193aSHong Zhang #undef __FUNCT__
5899123193aSHong Zhang #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIDense"
5909123193aSHong Zhang PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIDense(Mat A,Mat B,Mat C)
5919123193aSHong Zhang {
5929123193aSHong Zhang   PetscErrorCode ierr;
593f32d5d43SBarry Smith   Mat_MPIAIJ     *aij    = (Mat_MPIAIJ*)A->data;
594f32d5d43SBarry Smith   Mat_MPIDense   *bdense = (Mat_MPIDense*)B->data;
595f32d5d43SBarry Smith   Mat_MPIDense   *cdense = (Mat_MPIDense*)C->data;
596f32d5d43SBarry Smith   Mat            workB;
5979123193aSHong Zhang 
5989123193aSHong Zhang   PetscFunctionBegin;
599f32d5d43SBarry Smith   /* diagonal block of A times all local rows of B*/
600f32d5d43SBarry Smith   ierr = MatMatMultNumeric_SeqAIJ_SeqDense(aij->A,bdense->A,cdense->A);CHKERRQ(ierr);
601f32d5d43SBarry Smith 
602f32d5d43SBarry Smith   /* get off processor parts of B needed to complete the product */
6034324174eSBarry Smith   ierr = MatMPIDenseScatter(A,B,C,&workB);CHKERRQ(ierr);
604f32d5d43SBarry Smith 
605f32d5d43SBarry Smith   /* off-diagonal block of A times nonlocal rows of B */
606f32d5d43SBarry Smith   ierr = MatMatMultNumericAdd_SeqAIJ_SeqDense(aij->B,workB,cdense->A);CHKERRQ(ierr);
6079123193aSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
6089123193aSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
6099123193aSHong Zhang   PetscFunctionReturn(0);
6109123193aSHong Zhang }
611cf1a0672SHong Zhang 
6121634b032SHong Zhang #undef __FUNCT__
613f8efd71cSHong Zhang #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIAIJ"
614f8efd71cSHong Zhang PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C)
6151634b032SHong Zhang {
6161634b032SHong Zhang   PetscErrorCode ierr;
617cf1a0672SHong Zhang   Mat_MPIAIJ     *a   = (Mat_MPIAIJ*)A->data,*c=(Mat_MPIAIJ*)C->data;
618cf1a0672SHong Zhang   Mat_SeqAIJ     *ad  = (Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
619cf1a0672SHong Zhang   Mat_SeqAIJ     *cd  = (Mat_SeqAIJ*)(c->A)->data,*co=(Mat_SeqAIJ*)(c->B)->data;
620cf1a0672SHong Zhang   PetscInt       *adi = ad->i,*adj,*aoi=ao->i,*aoj;
621cf1a0672SHong Zhang   PetscScalar    *ada,*aoa,*cda=cd->a,*coa=co->a;
622cf1a0672SHong Zhang   Mat_SeqAIJ     *p_loc,*p_oth;
623cf1a0672SHong Zhang   PetscInt       *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pj;
62429825010SHong Zhang   PetscScalar    *pa_loc,*pa_oth,*pa,valtmp,*ca;
62529825010SHong Zhang   PetscInt       cm          = C->rmap->n,anz,pnz;
626cf1a0672SHong Zhang   Mat_PtAPMPI    *ptap       = c->ptap;
62729825010SHong Zhang   PetscScalar    *apa_sparse = ptap->apa;
628cf1a0672SHong Zhang   PetscInt       *api,*apj,*apJ,i,j,k,row;
62929825010SHong Zhang   PetscInt       cstart = C->cmap->rstart;
63029825010SHong Zhang   PetscInt       cdnz,conz,k0,k1,nextp;
631a7c7454dSHong Zhang   MPI_Comm       comm;
632a7c7454dSHong Zhang   PetscMPIInt    size;
6331634b032SHong Zhang 
6341634b032SHong Zhang   PetscFunctionBegin;
635a7c7454dSHong Zhang   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
636a7c7454dSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
637a7c7454dSHong Zhang 
638cf1a0672SHong Zhang   /* 1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
639cf1a0672SHong Zhang   /*-----------------------------------------------------*/
640cf1a0672SHong Zhang   /* update numerical values of P_oth and P_loc */
641b7f45c76SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
642cf1a0672SHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
643cf1a0672SHong Zhang 
644cf1a0672SHong Zhang   /* 2) compute numeric C_loc = A_loc*P = Ad*P_loc + Ao*P_oth */
645cf1a0672SHong Zhang   /*----------------------------------------------------------*/
646cf1a0672SHong Zhang   /* get data from symbolic products */
647cf1a0672SHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
648cf1a0672SHong Zhang   pi_loc = p_loc->i; pj_loc = p_loc->j; pa_loc = p_loc->a;
649a7c7454dSHong Zhang   if (size >1) {
650a7c7454dSHong Zhang     p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
651cf1a0672SHong Zhang     pi_oth = p_oth->i; pj_oth = p_oth->j; pa_oth = p_oth->a;
65220e3dc75SHong Zhang   } else {
65320e3dc75SHong Zhang     p_oth = NULL; pi_oth = NULL; pj_oth = NULL; pa_oth = NULL;
654a7c7454dSHong Zhang   }
655cf1a0672SHong Zhang 
65629825010SHong Zhang   api = ptap->api;
65729825010SHong Zhang   apj = ptap->apj;
658cf1a0672SHong Zhang   for (i=0; i<cm; i++) {
65929825010SHong Zhang     apJ = apj + api[i];
66029825010SHong Zhang 
661cf1a0672SHong Zhang     /* diagonal portion of A */
662cf1a0672SHong Zhang     anz = adi[i+1] - adi[i];
663cf1a0672SHong Zhang     adj = ad->j + adi[i];
664cf1a0672SHong Zhang     ada = ad->a + adi[i];
665cf1a0672SHong Zhang     for (j=0; j<anz; j++) {
666cf1a0672SHong Zhang       row = adj[j];
667cf1a0672SHong Zhang       pnz = pi_loc[row+1] - pi_loc[row];
668cf1a0672SHong Zhang       pj  = pj_loc + pi_loc[row];
669cf1a0672SHong Zhang       pa  = pa_loc + pi_loc[row];
67029825010SHong Zhang       /* perform sparse axpy */
671cf1a0672SHong Zhang       valtmp = ada[j];
67229825010SHong Zhang       nextp  = 0;
67329825010SHong Zhang       for (k=0; nextp<pnz; k++) {
67429825010SHong Zhang         if (apJ[k] == pj[nextp]) { /* column of AP == column of P */
67529825010SHong Zhang           apa_sparse[k] += valtmp*pa[nextp++];
67629825010SHong Zhang         }
6771634b032SHong Zhang       }
678cf1a0672SHong Zhang       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
679cf1a0672SHong Zhang     }
6801634b032SHong Zhang 
681cf1a0672SHong Zhang     /* off-diagonal portion of A */
682cf1a0672SHong Zhang     anz = aoi[i+1] - aoi[i];
683cf1a0672SHong Zhang     aoj = ao->j + aoi[i];
684cf1a0672SHong Zhang     aoa = ao->a + aoi[i];
685cf1a0672SHong Zhang     for (j=0; j<anz; j++) {
686cf1a0672SHong Zhang       row = aoj[j];
687cf1a0672SHong Zhang       pnz = pi_oth[row+1] - pi_oth[row];
688cf1a0672SHong Zhang       pj  = pj_oth + pi_oth[row];
689cf1a0672SHong Zhang       pa  = pa_oth + pi_oth[row];
69029825010SHong Zhang       /* perform sparse axpy */
691cf1a0672SHong Zhang       valtmp = aoa[j];
69229825010SHong Zhang       nextp  = 0;
69329825010SHong Zhang       for (k=0; nextp<pnz; k++) {
69429825010SHong Zhang         if (apJ[k] == pj[nextp]) { /* column of AP == column of P */
69529825010SHong Zhang           apa_sparse[k] += valtmp*pa[nextp++];
69629825010SHong Zhang         }
697cf1a0672SHong Zhang       }
698cf1a0672SHong Zhang       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
699cf1a0672SHong Zhang     }
7001634b032SHong Zhang 
701cf1a0672SHong Zhang     /* set values in C */
702cf1a0672SHong Zhang     cdnz = cd->i[i+1] - cd->i[i];
703cf1a0672SHong Zhang     conz = co->i[i+1] - co->i[i];
7041634b032SHong Zhang 
705cf1a0672SHong Zhang     /* 1st off-diagoanl part of C */
706cf1a0672SHong Zhang     ca = coa + co->i[i];
707cf1a0672SHong Zhang     k  = 0;
708cf1a0672SHong Zhang     for (k0=0; k0<conz; k0++) {
709cf1a0672SHong Zhang       if (apJ[k] >= cstart) break;
71029825010SHong Zhang       ca[k0]        = apa_sparse[k];
71129825010SHong Zhang       apa_sparse[k] = 0.0;
712cf1a0672SHong Zhang       k++;
713cf1a0672SHong Zhang     }
7141634b032SHong Zhang 
715cf1a0672SHong Zhang     /* diagonal part of C */
716cf1a0672SHong Zhang     ca = cda + cd->i[i];
717cf1a0672SHong Zhang     for (k1=0; k1<cdnz; k1++) {
71829825010SHong Zhang       ca[k1]        = apa_sparse[k];
71929825010SHong Zhang       apa_sparse[k] = 0.0;
720cf1a0672SHong Zhang       k++;
721cf1a0672SHong Zhang     }
722cf1a0672SHong Zhang 
723cf1a0672SHong Zhang     /* 2nd off-diagoanl part of C */
724cf1a0672SHong Zhang     ca = coa + co->i[i];
725cf1a0672SHong Zhang     for (; k0<conz; k0++) {
72629825010SHong Zhang       ca[k0]        = apa_sparse[k];
72729825010SHong Zhang       apa_sparse[k] = 0.0;
728cf1a0672SHong Zhang       k++;
729cf1a0672SHong Zhang     }
730cf1a0672SHong Zhang   }
731cf1a0672SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
732cf1a0672SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
733cf1a0672SHong Zhang   PetscFunctionReturn(0);
734cf1a0672SHong Zhang }
735cf1a0672SHong Zhang 
7360fc8cf34SHong Zhang /* same as MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(), except using LLCondensed to avoid O(BN) memory requirement */
737cf1a0672SHong Zhang #undef __FUNCT__
738b2405163SHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ"
739b2405163SHong Zhang PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
740cf1a0672SHong Zhang {
741cf1a0672SHong Zhang   PetscErrorCode     ierr;
742ce94432eSBarry Smith   MPI_Comm           comm;
743a7c7454dSHong Zhang   PetscMPIInt        size;
744cf1a0672SHong Zhang   Mat                Cmpi;
745cf1a0672SHong Zhang   Mat_PtAPMPI        *ptap;
7460298fd71SBarry Smith   PetscFreeSpaceList free_space = NULL,current_space=NULL;
747cf1a0672SHong Zhang   Mat_MPIAIJ         *a         = (Mat_MPIAIJ*)A->data,*c;
748cf1a0672SHong Zhang   Mat_SeqAIJ         *ad        = (Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_loc,*p_oth;
749cf1a0672SHong Zhang   PetscInt           *pi_loc,*pj_loc,*pi_oth,*pj_oth,*dnz,*onz;
750cf1a0672SHong Zhang   PetscInt           *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,rstart=A->rmap->rstart;
751f84c536bSHong Zhang   PetscInt           i,pnz,row,*api,*apj,*Jptr,apnz,nspacedouble=0,j,nzi,*lnk,apnz_max=0;
752cf1a0672SHong Zhang   PetscInt           am=A->rmap->n,pN=P->cmap->N,pn=P->cmap->n,pm=P->rmap->n;
753cf1a0672SHong Zhang   PetscInt           nlnk_max,armax,prmax;
754cf1a0672SHong Zhang   PetscReal          afill;
755cf1a0672SHong Zhang   PetscScalar        *apa;
756cf1a0672SHong Zhang 
757cf1a0672SHong Zhang   PetscFunctionBegin;
758ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
759a7c7454dSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
760a7c7454dSHong Zhang 
761cf1a0672SHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
762b00a9115SJed Brown   ierr = PetscNew(&ptap);CHKERRQ(ierr);
763cf1a0672SHong Zhang 
764cf1a0672SHong Zhang   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
765b7f45c76SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
76619f0e57aSHong Zhang 
767cf1a0672SHong Zhang   /* get P_loc by taking all local rows of P */
768cf1a0672SHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
769cf1a0672SHong Zhang 
770cf1a0672SHong Zhang   p_loc  = (Mat_SeqAIJ*)(ptap->P_loc)->data;
771cf1a0672SHong Zhang   pi_loc = p_loc->i; pj_loc = p_loc->j;
772a7c7454dSHong Zhang   if (size > 1) {
773a7c7454dSHong Zhang     p_oth  = (Mat_SeqAIJ*)(ptap->P_oth)->data;
77420e3dc75SHong Zhang     pi_oth = p_oth->i; pj_oth = p_oth->j;
775968382fdSHong Zhang   } else {
77620e3dc75SHong Zhang     p_oth  = NULL;
77720e3dc75SHong Zhang     pi_oth = NULL; pj_oth = NULL;
778a7c7454dSHong Zhang   }
779cf1a0672SHong Zhang 
780cf1a0672SHong Zhang   /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */
781cf1a0672SHong Zhang   /*-------------------------------------------------------------------*/
782854ce69bSBarry Smith   ierr      = PetscMalloc1(am+2,&api);CHKERRQ(ierr);
783cf1a0672SHong Zhang   ptap->api = api;
784cf1a0672SHong Zhang   api[0]    = 0;
785cf1a0672SHong Zhang 
786cf1a0672SHong Zhang   /* create and initialize a linked list */
787cf1a0672SHong Zhang   armax = ad->rmax+ao->rmax;
788a7c7454dSHong Zhang   if (size >1) {
789cf1a0672SHong Zhang     prmax = PetscMax(p_loc->rmax,p_oth->rmax);
790a7c7454dSHong Zhang   } else {
791a7c7454dSHong Zhang     prmax = p_loc->rmax;
792a7c7454dSHong Zhang   }
793cf1a0672SHong Zhang   nlnk_max = armax*prmax;
794cf1a0672SHong Zhang   if (!nlnk_max || nlnk_max > pN) nlnk_max = pN;
795f84c536bSHong Zhang   ierr = PetscLLCondensedCreate_Scalable(nlnk_max,&lnk);CHKERRQ(ierr);
796cf1a0672SHong Zhang 
797cf1a0672SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */
798cf1a0672SHong Zhang   ierr = PetscFreeSpaceGet((PetscInt)(fill*(adi[am]+aoi[am]+pi_loc[pm])),&free_space);CHKERRQ(ierr);
7992205254eSKarl Rupp 
800cf1a0672SHong Zhang   current_space = free_space;
801cf1a0672SHong Zhang 
802cf1a0672SHong Zhang   ierr = MatPreallocateInitialize(comm,am,pn,dnz,onz);CHKERRQ(ierr);
803cf1a0672SHong Zhang   for (i=0; i<am; i++) {
804cf1a0672SHong Zhang     /* diagonal portion of A */
805cf1a0672SHong Zhang     nzi = adi[i+1] - adi[i];
806cf1a0672SHong Zhang     for (j=0; j<nzi; j++) {
807cf1a0672SHong Zhang       row  = *adj++;
808cf1a0672SHong Zhang       pnz  = pi_loc[row+1] - pi_loc[row];
809cf1a0672SHong Zhang       Jptr = pj_loc + pi_loc[row];
810cf1a0672SHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
811f84c536bSHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr);
812cf1a0672SHong Zhang     }
813cf1a0672SHong Zhang     /* off-diagonal portion of A */
814cf1a0672SHong Zhang     nzi = aoi[i+1] - aoi[i];
815cf1a0672SHong Zhang     for (j=0; j<nzi; j++) {
816cf1a0672SHong Zhang       row  = *aoj++;
817cf1a0672SHong Zhang       pnz  = pi_oth[row+1] - pi_oth[row];
818cf1a0672SHong Zhang       Jptr = pj_oth + pi_oth[row];
819f84c536bSHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr);
820cf1a0672SHong Zhang     }
821cf1a0672SHong Zhang 
822f84c536bSHong Zhang     apnz     = *lnk;
823cf1a0672SHong Zhang     api[i+1] = api[i] + apnz;
824cf1a0672SHong Zhang     if (apnz > apnz_max) apnz_max = apnz;
825cf1a0672SHong Zhang 
826cf1a0672SHong Zhang     /* if free space is not available, double the total space in the list */
827cf1a0672SHong Zhang     if (current_space->local_remaining<apnz) {
828cf1a0672SHong Zhang       ierr = PetscFreeSpaceGet(apnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
829cf1a0672SHong Zhang       nspacedouble++;
830cf1a0672SHong Zhang     }
831cf1a0672SHong Zhang 
832cf1a0672SHong Zhang     /* Copy data into free space, then initialize lnk */
833f84c536bSHong Zhang     ierr = PetscLLCondensedClean_Scalable(apnz,current_space->array,lnk);CHKERRQ(ierr);
834cf1a0672SHong Zhang     ierr = MatPreallocateSet(i+rstart,apnz,current_space->array,dnz,onz);CHKERRQ(ierr);
8352205254eSKarl Rupp 
836cf1a0672SHong Zhang     current_space->array           += apnz;
837cf1a0672SHong Zhang     current_space->local_used      += apnz;
838cf1a0672SHong Zhang     current_space->local_remaining -= apnz;
839cf1a0672SHong Zhang   }
840cf1a0672SHong Zhang 
841cf1a0672SHong Zhang   /* Allocate space for apj, initialize apj, and */
842cf1a0672SHong Zhang   /* destroy list of free space and other temporary array(s) */
843854ce69bSBarry Smith   ierr = PetscMalloc1(api[am]+1,&ptap->apj);CHKERRQ(ierr);
844cf1a0672SHong Zhang   apj  = ptap->apj;
845cf1a0672SHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,ptap->apj);CHKERRQ(ierr);
846f84c536bSHong Zhang   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr);
847cf1a0672SHong Zhang 
848cf1a0672SHong Zhang   /* create and assemble symbolic parallel matrix Cmpi */
849cf1a0672SHong Zhang   /*----------------------------------------------------*/
850cf1a0672SHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
851cf1a0672SHong Zhang   ierr = MatSetSizes(Cmpi,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
85233d57670SJed Brown   ierr = MatSetBlockSizesFromMats(Cmpi,A,P);CHKERRQ(ierr);
853cf1a0672SHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
854cf1a0672SHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
855cf1a0672SHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
856cf1a0672SHong Zhang 
85729825010SHong Zhang   /* malloc apa for assembly Cmpi */
8581795a4d1SJed Brown   ierr = PetscCalloc1(apnz_max,&apa);CHKERRQ(ierr);
8592205254eSKarl Rupp 
86029825010SHong Zhang   ptap->apa = apa;
861cf1a0672SHong Zhang   for (i=0; i<am; i++) {
862cf1a0672SHong Zhang     row  = i + rstart;
863cf1a0672SHong Zhang     apnz = api[i+1] - api[i];
864cf1a0672SHong Zhang     ierr = MatSetValues(Cmpi,1,&row,apnz,apj,apa,INSERT_VALUES);CHKERRQ(ierr);
865cf1a0672SHong Zhang     apj += apnz;
866cf1a0672SHong Zhang   }
867cf1a0672SHong Zhang   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
868cf1a0672SHong Zhang   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
869cf1a0672SHong Zhang 
870cf1a0672SHong Zhang   ptap->destroy             = Cmpi->ops->destroy;
871cf1a0672SHong Zhang   ptap->duplicate           = Cmpi->ops->duplicate;
872f8efd71cSHong Zhang   Cmpi->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ;
873cf1a0672SHong Zhang   Cmpi->ops->destroy        = MatDestroy_MPIAIJ_MatMatMult;
874cf1a0672SHong Zhang   Cmpi->ops->duplicate      = MatDuplicate_MPIAIJ_MatMatMult;
875cf1a0672SHong Zhang 
876cf1a0672SHong Zhang   /* attach the supporting struct to Cmpi for reuse */
877cf1a0672SHong Zhang   c       = (Mat_MPIAIJ*)Cmpi->data;
878cf1a0672SHong Zhang   c->ptap = ptap;
879cf1a0672SHong Zhang 
880cf1a0672SHong Zhang   *C = Cmpi;
881cf1a0672SHong Zhang 
882cf1a0672SHong Zhang   /* set MatInfo */
883118727c9SMark F. Adams   afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1) + 1.e-5;
884cf1a0672SHong Zhang   if (afill < 1.0) afill = 1.0;
885cf1a0672SHong Zhang   Cmpi->info.mallocs           = nspacedouble;
886cf1a0672SHong Zhang   Cmpi->info.fill_ratio_given  = fill;
887cf1a0672SHong Zhang   Cmpi->info.fill_ratio_needed = afill;
888cf1a0672SHong Zhang 
889cf1a0672SHong Zhang #if defined(PETSC_USE_INFO)
890cf1a0672SHong Zhang   if (api[am]) {
89157622a8eSBarry Smith     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
89257622a8eSBarry Smith     ierr = PetscInfo1(Cmpi,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr);
893cf1a0672SHong Zhang   } else {
894cf1a0672SHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
895cf1a0672SHong Zhang   }
896cf1a0672SHong Zhang #endif
8971634b032SHong Zhang   PetscFunctionReturn(0);
8981634b032SHong Zhang }
899f96d37fcSHong Zhang 
900f96d37fcSHong Zhang /*-------------------------------------------------------------------------*/
901f96d37fcSHong Zhang #undef __FUNCT__
902f96d37fcSHong Zhang #define __FUNCT__ "MatTransposeMatMult_MPIAIJ_MPIAIJ"
903c5af039cSHong Zhang PetscErrorCode MatTransposeMatMult_MPIAIJ_MPIAIJ(Mat P,Mat A,MatReuse scall,PetscReal fill,Mat *C)
904f96d37fcSHong Zhang {
905f96d37fcSHong Zhang   PetscErrorCode ierr;
906c216dbf3SHong Zhang   const char     *algTypes[3] = {"scalable","nonscalable","matmatmult"};
907c216dbf3SHong Zhang   PetscInt       alg=0; /* set default algorithm */
908f96d37fcSHong Zhang 
909f96d37fcSHong Zhang   PetscFunctionBegin;
910c216dbf3SHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
911c216dbf3SHong Zhang     ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr);
912c216dbf3SHong Zhang     ierr = PetscOptionsEList("-mattransposematmult_via","Algorithmic approach","MatTransposeMatMult",algTypes,3,algTypes[0],&alg,NULL);CHKERRQ(ierr);
913c216dbf3SHong Zhang     ierr = PetscOptionsEnd();CHKERRQ(ierr);
914c216dbf3SHong Zhang 
915c216dbf3SHong Zhang     ierr = PetscLogEventBegin(MAT_TransposeMatMultSymbolic,P,A,0,0);CHKERRQ(ierr);
916c216dbf3SHong Zhang     switch (alg) {
917c216dbf3SHong Zhang     case 1:
9182bbb1c24SHong Zhang       ierr = MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(P,A,fill,C);CHKERRQ(ierr);
919c216dbf3SHong Zhang       break;
920c216dbf3SHong Zhang     case 2:
921275476c6SMatthew G. Knepley     {
922ab1b013aSHong Zhang       Mat         Pt;
923ab1b013aSHong Zhang       Mat_PtAPMPI *ptap;
924ab1b013aSHong Zhang       Mat_MPIAIJ  *c;
925ab1b013aSHong Zhang       ierr = MatTranspose(P,MAT_INITIAL_MATRIX,&Pt);CHKERRQ(ierr);
926ab1b013aSHong Zhang       ierr = MatMatMult(Pt,A,MAT_INITIAL_MATRIX,fill,C);CHKERRQ(ierr);
927ab1b013aSHong Zhang       c        = (Mat_MPIAIJ*)(*C)->data;
928ab1b013aSHong Zhang       ptap     = c->ptap;
929ab1b013aSHong Zhang       ptap->Pt = Pt;
930c216dbf3SHong Zhang       (*C)->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult;
931c216dbf3SHong Zhang       PetscFunctionReturn(0);
932275476c6SMatthew G. Knepley     }
933c216dbf3SHong Zhang       break;
934c216dbf3SHong Zhang     default:
9356da69ca6SHong Zhang       ierr = MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(P,A,fill,C);CHKERRQ(ierr);
936c216dbf3SHong Zhang       break;
937c216dbf3SHong Zhang     }
938c216dbf3SHong Zhang     ierr = PetscLogEventEnd(MAT_TransposeMatMultSymbolic,P,A,0,0);CHKERRQ(ierr);
939c216dbf3SHong Zhang   }
940c216dbf3SHong Zhang   ierr = PetscLogEventBegin(MAT_TransposeMatMultNumeric,P,A,0,0);CHKERRQ(ierr);
941c216dbf3SHong Zhang   ierr = (*(*C)->ops->mattransposemultnumeric)(P,A,*C);CHKERRQ(ierr);
942c216dbf3SHong Zhang   ierr = PetscLogEventEnd(MAT_TransposeMatMultNumeric,P,A,0,0);CHKERRQ(ierr);
943ab1b013aSHong Zhang   PetscFunctionReturn(0);
944ab1b013aSHong Zhang }
945ab1b013aSHong Zhang 
946c216dbf3SHong Zhang /* This routine only works when scall=MAT_REUSE_MATRIX! */
947c216dbf3SHong Zhang #undef __FUNCT__
948c216dbf3SHong Zhang #define __FUNCT__ "MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult"
949c216dbf3SHong Zhang PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult(Mat P,Mat A,Mat C)
950c216dbf3SHong Zhang {
951c216dbf3SHong Zhang   PetscErrorCode ierr;
9522bbb1c24SHong Zhang   Mat_MPIAIJ     *c=(Mat_MPIAIJ*)C->data;
9532bbb1c24SHong Zhang   Mat_PtAPMPI    *ptap= c->ptap;
9542bbb1c24SHong Zhang   Mat            Pt=ptap->Pt;
955c216dbf3SHong Zhang 
956c216dbf3SHong Zhang   PetscFunctionBegin;
957c216dbf3SHong Zhang   ierr = MatTranspose(P,MAT_REUSE_MATRIX,&Pt);CHKERRQ(ierr);
958c216dbf3SHong Zhang   ierr = MatMatMultNumeric(Pt,A,C);CHKERRQ(ierr);
959f96d37fcSHong Zhang   PetscFunctionReturn(0);
960f96d37fcSHong Zhang }
961f96d37fcSHong Zhang 
9626da69ca6SHong Zhang /* Non-scalable version, use dense axpy */
963f96d37fcSHong Zhang #undef __FUNCT__
9646da69ca6SHong Zhang #define __FUNCT__ "MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable"
9656da69ca6SHong Zhang PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable(Mat P,Mat A,Mat C)
966f96d37fcSHong Zhang {
967c5af039cSHong Zhang   PetscErrorCode      ierr;
968c5af039cSHong Zhang   Mat_Merge_SeqsToMPI *merge;
969497f5370SHong Zhang   Mat_MPIAIJ          *p =(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
970c5af039cSHong Zhang   Mat_SeqAIJ          *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data;
971c5af039cSHong Zhang   Mat_PtAPMPI         *ptap;
972d6ab1dc0SHong Zhang   PetscInt            *adj,*aJ;
973497f5370SHong Zhang   PetscInt            i,j,k,anz,pnz,row,*cj;
974e745005bSHong Zhang   MatScalar           *ada,*aval,*ca,valtmp;
975c5af039cSHong Zhang   PetscInt            am  =A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n;
976ce94432eSBarry Smith   MPI_Comm            comm;
977c5af039cSHong Zhang   PetscMPIInt         size,rank,taga,*len_s;
978c5af039cSHong Zhang   PetscInt            *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci;
979c5af039cSHong Zhang   PetscInt            **buf_ri,**buf_rj;
980c5af039cSHong Zhang   PetscInt            cnz=0,*bj_i,*bi,*bj,bnz,nextcj;  /* bi,bj,ba: local array of C(mpi mat) */
981c5af039cSHong Zhang   MPI_Request         *s_waits,*r_waits;
982c5af039cSHong Zhang   MPI_Status          *status;
983c5af039cSHong Zhang   MatScalar           **abuf_r,*ba_i,*pA,*coa,*ba;
984d6ab1dc0SHong Zhang   PetscInt            *ai,*aj,*coi,*coj;
985a2ea699eSBarry Smith   PetscInt            *poJ,*pdJ;
986c5af039cSHong Zhang   Mat                 A_loc;
987c5af039cSHong Zhang   Mat_SeqAIJ          *a_loc;
988f96d37fcSHong Zhang 
989f96d37fcSHong Zhang   PetscFunctionBegin;
990ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr);
991c5af039cSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
992c5af039cSHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
993c5af039cSHong Zhang 
994c5af039cSHong Zhang   ptap  = c->ptap;
995c5af039cSHong Zhang   merge = ptap->merge;
996c5af039cSHong Zhang 
997c5af039cSHong Zhang   /* 2) compute numeric C_seq = P_loc^T*A_loc*P - dominating part */
998c5af039cSHong Zhang   /*--------------------------------------------------------------*/
999c5af039cSHong Zhang   /* get data from symbolic products */
1000c5af039cSHong Zhang   coi  = merge->coi; coj = merge->coj;
1001854ce69bSBarry Smith   ierr = PetscCalloc1(coi[pon]+1,&coa);CHKERRQ(ierr);
1002c5af039cSHong Zhang 
1003c5af039cSHong Zhang   bi     = merge->bi; bj = merge->bj;
1004c5af039cSHong Zhang   owners = merge->rowmap->range;
1005854ce69bSBarry Smith   ierr   = PetscCalloc1(bi[cm]+1,&ba);CHKERRQ(ierr);
1006c5af039cSHong Zhang 
1007c5af039cSHong Zhang   /* get A_loc by taking all local rows of A */
1008c5af039cSHong Zhang   A_loc = ptap->A_loc;
1009c5af039cSHong Zhang   ierr  = MatMPIAIJGetLocalMat(A,MAT_REUSE_MATRIX,&A_loc);CHKERRQ(ierr);
1010c5af039cSHong Zhang   a_loc = (Mat_SeqAIJ*)(A_loc)->data;
1011d6ab1dc0SHong Zhang   ai    = a_loc->i;
1012d6ab1dc0SHong Zhang   aj    = a_loc->j;
1013c5af039cSHong Zhang 
1014854ce69bSBarry Smith   ierr = PetscCalloc1(A->cmap->N,&aval);CHKERRQ(ierr); /* non-scalable!!! */
1015c5af039cSHong Zhang 
1016c5af039cSHong Zhang   for (i=0; i<am; i++) {
1017e745005bSHong Zhang     /* 2-a) put A[i,:] to dense array aval */
1018d6ab1dc0SHong Zhang     anz = ai[i+1] - ai[i];
1019d6ab1dc0SHong Zhang     adj = aj + ai[i];
1020d6ab1dc0SHong Zhang     ada = a_loc->a + ai[i];
1021c5af039cSHong Zhang     for (j=0; j<anz; j++) {
1022e745005bSHong Zhang       aval[adj[j]] = ada[j];
1023c5af039cSHong Zhang     }
1024c5af039cSHong Zhang 
1025c5af039cSHong Zhang     /* 2-b) Compute Cseq = P_loc[i,:]^T*A[i,:] using outer product */
1026c5af039cSHong Zhang     /*--------------------------------------------------------------*/
1027c5af039cSHong Zhang     /* put the value into Co=(p->B)^T*A (off-diagonal part, send to others) */
1028c5af039cSHong Zhang     pnz = po->i[i+1] - po->i[i];
1029c5af039cSHong Zhang     poJ = po->j + po->i[i];
1030c5af039cSHong Zhang     pA  = po->a + po->i[i];
1031c5af039cSHong Zhang     for (j=0; j<pnz; j++) {
1032c5af039cSHong Zhang       row = poJ[j];
1033c5af039cSHong Zhang       cnz = coi[row+1] - coi[row];
1034c5af039cSHong Zhang       cj  = coj + coi[row];
1035c5af039cSHong Zhang       ca  = coa + coi[row];
1036c5af039cSHong Zhang       /* perform dense axpy */
1037c5af039cSHong Zhang       valtmp = pA[j];
1038c5af039cSHong Zhang       for (k=0; k<cnz; k++) {
1039e745005bSHong Zhang         ca[k] += valtmp*aval[cj[k]];
1040c5af039cSHong Zhang       }
1041c5af039cSHong Zhang       ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1042c5af039cSHong Zhang     }
1043c5af039cSHong Zhang 
1044c5af039cSHong Zhang     /* put the value into Cd (diagonal part) */
1045c5af039cSHong Zhang     pnz = pd->i[i+1] - pd->i[i];
1046c5af039cSHong Zhang     pdJ = pd->j + pd->i[i];
1047c5af039cSHong Zhang     pA  = pd->a + pd->i[i];
1048c5af039cSHong Zhang     for (j=0; j<pnz; j++) {
1049c5af039cSHong Zhang       row = pdJ[j];
1050c5af039cSHong Zhang       cnz = bi[row+1] - bi[row];
1051c5af039cSHong Zhang       cj  = bj + bi[row];
1052c5af039cSHong Zhang       ca  = ba + bi[row];
1053c5af039cSHong Zhang       /* perform dense axpy */
1054c5af039cSHong Zhang       valtmp = pA[j];
1055c5af039cSHong Zhang       for (k=0; k<cnz; k++) {
1056e745005bSHong Zhang         ca[k] += valtmp*aval[cj[k]];
1057c5af039cSHong Zhang       }
1058c5af039cSHong Zhang       ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1059c5af039cSHong Zhang     }
1060c5af039cSHong Zhang 
1061d6ab1dc0SHong Zhang     /* zero the current row of Pt*A */
1062d6ab1dc0SHong Zhang     aJ = aj + ai[i];
1063e745005bSHong Zhang     for (k=0; k<anz; k++) aval[aJ[k]] = 0.0;
1064c5af039cSHong Zhang   }
1065c5af039cSHong Zhang 
1066c5af039cSHong Zhang   /* 3) send and recv matrix values coa */
1067c5af039cSHong Zhang   /*------------------------------------*/
1068c5af039cSHong Zhang   buf_ri = merge->buf_ri;
1069c5af039cSHong Zhang   buf_rj = merge->buf_rj;
1070c5af039cSHong Zhang   len_s  = merge->len_s;
1071c5af039cSHong Zhang   ierr   = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr);
1072c5af039cSHong Zhang   ierr   = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr);
1073c5af039cSHong Zhang 
1074dcca6d9dSJed Brown   ierr = PetscMalloc2(merge->nsend+1,&s_waits,size,&status);CHKERRQ(ierr);
1075c5af039cSHong Zhang   for (proc=0,k=0; proc<size; proc++) {
1076c5af039cSHong Zhang     if (!len_s[proc]) continue;
1077c5af039cSHong Zhang     i    = merge->owners_co[proc];
1078c5af039cSHong Zhang     ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr);
1079c5af039cSHong Zhang     k++;
1080c5af039cSHong Zhang   }
1081c5af039cSHong Zhang   if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);}
1082c5af039cSHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);}
1083c5af039cSHong Zhang 
1084c5af039cSHong Zhang   ierr = PetscFree2(s_waits,status);CHKERRQ(ierr);
1085c5af039cSHong Zhang   ierr = PetscFree(r_waits);CHKERRQ(ierr);
1086c5af039cSHong Zhang   ierr = PetscFree(coa);CHKERRQ(ierr);
1087c5af039cSHong Zhang 
1088c5af039cSHong Zhang   /* 4) insert local Cseq and received values into Cmpi */
1089c5af039cSHong Zhang   /*----------------------------------------------------*/
1090dcca6d9dSJed Brown   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
1091c5af039cSHong Zhang   for (k=0; k<merge->nrecv; k++) {
1092c36aecf5SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1093c5af039cSHong Zhang     nrows       = *(buf_ri_k[k]);
1094c5af039cSHong Zhang     nextrow[k]  = buf_ri_k[k]+1;  /* next row number of k-th recved i-structure */
1095c5af039cSHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
1096c5af039cSHong Zhang   }
1097c5af039cSHong Zhang 
1098c5af039cSHong Zhang   for (i=0; i<cm; i++) {
1099c5af039cSHong Zhang     row  = owners[rank] + i; /* global row index of C_seq */
1100c5af039cSHong Zhang     bj_i = bj + bi[i];  /* col indices of the i-th row of C */
1101c5af039cSHong Zhang     ba_i = ba + bi[i];
1102c5af039cSHong Zhang     bnz  = bi[i+1] - bi[i];
1103c5af039cSHong Zhang     /* add received vals into ba */
1104c5af039cSHong Zhang     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1105c5af039cSHong Zhang       /* i-th row */
1106c5af039cSHong Zhang       if (i == *nextrow[k]) {
1107c5af039cSHong Zhang         cnz    = *(nextci[k]+1) - *nextci[k];
1108c5af039cSHong Zhang         cj     = buf_rj[k] + *(nextci[k]);
1109c5af039cSHong Zhang         ca     = abuf_r[k] + *(nextci[k]);
1110c5af039cSHong Zhang         nextcj = 0;
1111c5af039cSHong Zhang         for (j=0; nextcj<cnz; j++) {
1112c5af039cSHong Zhang           if (bj_i[j] == cj[nextcj]) { /* bcol == ccol */
1113c5af039cSHong Zhang             ba_i[j] += ca[nextcj++];
1114c5af039cSHong Zhang           }
1115c5af039cSHong Zhang         }
1116c5af039cSHong Zhang         nextrow[k]++; nextci[k]++;
1117c5af039cSHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1118c5af039cSHong Zhang       }
1119c5af039cSHong Zhang     }
1120c5af039cSHong Zhang     ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr);
1121c5af039cSHong Zhang   }
1122c5af039cSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1123c5af039cSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1124c5af039cSHong Zhang 
1125c5af039cSHong Zhang   ierr = PetscFree(ba);CHKERRQ(ierr);
1126c5af039cSHong Zhang   ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr);
1127c5af039cSHong Zhang   ierr = PetscFree(abuf_r);CHKERRQ(ierr);
1128c5af039cSHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1129e745005bSHong Zhang   ierr = PetscFree(aval);CHKERRQ(ierr);
1130f96d37fcSHong Zhang   PetscFunctionReturn(0);
1131f96d37fcSHong Zhang }
1132f96d37fcSHong Zhang 
1133c9ba551fSBarry Smith PetscErrorCode MatDuplicate_MPIAIJ_MatPtAP(Mat, MatDuplicateOption,Mat*);
1134c5af039cSHong Zhang /* This routine is modified from MatPtAPSymbolic_MPIAIJ_MPIAIJ() */
1135f96d37fcSHong Zhang #undef __FUNCT__
11362bbb1c24SHong Zhang #define __FUNCT__ "MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable"
11372bbb1c24SHong Zhang PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat P,Mat A,PetscReal fill,Mat *C)
1138f96d37fcSHong Zhang {
1139f96d37fcSHong Zhang   PetscErrorCode      ierr;
11404e938277SHong Zhang   Mat                 Cmpi,A_loc,POt,PDt;
1141f96d37fcSHong Zhang   Mat_PtAPMPI         *ptap;
11420298fd71SBarry Smith   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
1143f96d37fcSHong Zhang   Mat_MPIAIJ          *p        =(Mat_MPIAIJ*)P->data,*c;
1144f96d37fcSHong Zhang   PetscInt            *pdti,*pdtj,*poti,*potj,*ptJ;
1145f96d37fcSHong Zhang   PetscInt            nnz;
11464e938277SHong Zhang   PetscInt            *lnk,*owners_co,*coi,*coj,i,k,pnz,row;
1147497f5370SHong Zhang   PetscInt            am=A->rmap->n,pn=P->cmap->n;
1148f96d37fcSHong Zhang   PetscBT             lnkbt;
1149ce94432eSBarry Smith   MPI_Comm            comm;
1150f96d37fcSHong Zhang   PetscMPIInt         size,rank,tagi,tagj,*len_si,*len_s,*len_ri;
1151f96d37fcSHong Zhang   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
1152f96d37fcSHong Zhang   PetscInt            len,proc,*dnz,*onz,*owners;
1153f96d37fcSHong Zhang   PetscInt            nzi,*bi,*bj;
1154f96d37fcSHong Zhang   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
1155f96d37fcSHong Zhang   MPI_Request         *swaits,*rwaits;
1156f96d37fcSHong Zhang   MPI_Status          *sstatus,rstatus;
1157f96d37fcSHong Zhang   Mat_Merge_SeqsToMPI *merge;
1158d6ab1dc0SHong Zhang   PetscInt            *ai,*aj,*Jptr,anz,*prmap=p->garray,pon,nspacedouble=0,j;
1159f96d37fcSHong Zhang   PetscReal           afill  =1.0,afill_tmp;
1160c36aecf5SHong Zhang   PetscInt            rstart = P->cmap->rstart,rmax,aN=A->cmap->N,Crmax;
1161f96d37fcSHong Zhang   PetscScalar         *vals;
11624e938277SHong Zhang   Mat_SeqAIJ          *a_loc, *pdt,*pot;
1163f96d37fcSHong Zhang 
1164f96d37fcSHong Zhang   PetscFunctionBegin;
1165ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
1166c5af039cSHong Zhang   /* check if matrix local sizes are compatible */
1167c5af039cSHong Zhang   if (A->rmap->rstart != P->rmap->rstart || A->rmap->rend != P->rmap->rend) {
1168c5af039cSHong Zhang     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);
1169c5af039cSHong Zhang   }
1170c5af039cSHong Zhang 
1171f96d37fcSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1172f96d37fcSHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1173f96d37fcSHong Zhang 
1174f96d37fcSHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
1175b00a9115SJed Brown   ierr = PetscNew(&ptap);CHKERRQ(ierr);
1176f96d37fcSHong Zhang 
1177f96d37fcSHong Zhang   /* get A_loc by taking all local rows of A */
1178f96d37fcSHong Zhang   ierr = MatMPIAIJGetLocalMat(A,MAT_INITIAL_MATRIX,&A_loc);CHKERRQ(ierr);
11792205254eSKarl Rupp 
1180c5af039cSHong Zhang   ptap->A_loc = A_loc;
11812205254eSKarl Rupp 
11821c7d5954SHong Zhang   a_loc = (Mat_SeqAIJ*)(A_loc)->data;
1183d6ab1dc0SHong Zhang   ai    = a_loc->i;
1184d6ab1dc0SHong Zhang   aj    = a_loc->j;
1185f96d37fcSHong Zhang 
1186f96d37fcSHong Zhang   /* determine symbolic Co=(p->B)^T*A - send to others */
1187f96d37fcSHong Zhang   /*----------------------------------------------------*/
11884e938277SHong Zhang   ierr = MatTransposeSymbolic_SeqAIJ(p->A,&PDt);CHKERRQ(ierr);
11894e938277SHong Zhang   pdt  = (Mat_SeqAIJ*)PDt->data;
11904e938277SHong Zhang   pdti = pdt->i; pdtj = pdt->j;
11914e938277SHong Zhang 
11924e938277SHong Zhang   ierr = MatTransposeSymbolic_SeqAIJ(p->B,&POt);CHKERRQ(ierr);
11934e938277SHong Zhang   pot  = (Mat_SeqAIJ*)POt->data;
11944e938277SHong Zhang   poti = pot->i; potj = pot->j;
1195f96d37fcSHong Zhang 
1196f96d37fcSHong Zhang   /* then, compute symbolic Co = (p->B)^T*A */
11972205254eSKarl Rupp   pon    = (p->B)->cmap->n; /* total num of rows to be sent to other processors >= (num of nonzero rows of C_seq) - pn */
1198854ce69bSBarry Smith   ierr   = PetscMalloc1(pon+1,&coi);CHKERRQ(ierr);
1199f96d37fcSHong Zhang   coi[0] = 0;
1200f96d37fcSHong Zhang 
1201f96d37fcSHong Zhang   /* set initial free space to be fill*(nnz(p->B) + nnz(A)) */
1202d6ab1dc0SHong Zhang   nnz           = fill*(poti[pon] + ai[am]);
1203a2ea699eSBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
1204f96d37fcSHong Zhang   current_space = free_space;
120519f0e57aSHong Zhang 
1206f96d37fcSHong Zhang   /* create and initialize a linked list */
12074e938277SHong Zhang   i     = PetscMax(pdt->rmax,pot->rmax);
1208c36aecf5SHong Zhang   Crmax = i*a_loc->rmax*size;
12094e938277SHong Zhang   if (!Crmax || Crmax > aN) Crmax = aN;
12104e938277SHong Zhang   ierr = PetscLLCondensedCreate(Crmax,aN,&lnk,&lnkbt);CHKERRQ(ierr);
1211f96d37fcSHong Zhang 
1212f96d37fcSHong Zhang   for (i=0; i<pon; i++) {
1213f96d37fcSHong Zhang     pnz = poti[i+1] - poti[i];
1214f96d37fcSHong Zhang     ptJ = potj + poti[i];
1215f96d37fcSHong Zhang     for (j=0; j<pnz; j++) {
1216f96d37fcSHong Zhang       row  = ptJ[j]; /* row of A_loc == col of Pot */
1217d6ab1dc0SHong Zhang       anz  = ai[row+1] - ai[row];
1218d6ab1dc0SHong Zhang       Jptr = aj + ai[row];
1219f96d37fcSHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
1220d6ab1dc0SHong Zhang       ierr = PetscLLCondensedAddSorted(anz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1221f96d37fcSHong Zhang     }
12224e938277SHong Zhang     nnz = lnk[0];
1223f96d37fcSHong Zhang 
1224f96d37fcSHong Zhang     /* If free space is not available, double the total space in the list */
1225f96d37fcSHong Zhang     if (current_space->local_remaining<nnz) {
1226f96d37fcSHong Zhang       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
1227f96d37fcSHong Zhang       nspacedouble++;
1228f96d37fcSHong Zhang     }
1229f96d37fcSHong Zhang 
1230f96d37fcSHong Zhang     /* Copy data into free space, and zero out denserows */
12314e938277SHong Zhang     ierr = PetscLLCondensedClean(aN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
12322205254eSKarl Rupp 
1233f96d37fcSHong Zhang     current_space->array           += nnz;
1234f96d37fcSHong Zhang     current_space->local_used      += nnz;
1235f96d37fcSHong Zhang     current_space->local_remaining -= nnz;
12362205254eSKarl Rupp 
1237f96d37fcSHong Zhang     coi[i+1] = coi[i] + nnz;
1238f96d37fcSHong Zhang   }
1239f96d37fcSHong Zhang 
1240854ce69bSBarry Smith   ierr = PetscMalloc1(coi[pon]+1,&coj);CHKERRQ(ierr);
1241f96d37fcSHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr);
12422205254eSKarl Rupp 
1243118727c9SMark F. Adams   afill_tmp = (PetscReal)coi[pon]/(poti[pon] + ai[am]+1);
1244f96d37fcSHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
1245f96d37fcSHong Zhang 
1246f96d37fcSHong Zhang   /* send j-array (coj) of Co to other processors */
1247f96d37fcSHong Zhang   /*----------------------------------------------*/
1248f96d37fcSHong Zhang   /* determine row ownership */
1249b00a9115SJed Brown   ierr = PetscNew(&merge);CHKERRQ(ierr);
1250f96d37fcSHong Zhang   ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr);
12512205254eSKarl Rupp 
1252f96d37fcSHong Zhang   merge->rowmap->n  = pn;
1253f96d37fcSHong Zhang   merge->rowmap->bs = 1;
12542205254eSKarl Rupp 
1255f96d37fcSHong Zhang   ierr   = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr);
1256f96d37fcSHong Zhang   owners = merge->rowmap->range;
1257f96d37fcSHong Zhang 
1258f96d37fcSHong Zhang   /* determine the number of messages to send, their lengths */
12591795a4d1SJed Brown   ierr = PetscCalloc1(size,&len_si);CHKERRQ(ierr);
1260785e854fSJed Brown   ierr = PetscMalloc1(size,&merge->len_s);CHKERRQ(ierr);
12612205254eSKarl Rupp 
1262f96d37fcSHong Zhang   len_s        = merge->len_s;
1263f96d37fcSHong Zhang   merge->nsend = 0;
1264f96d37fcSHong Zhang 
1265854ce69bSBarry Smith   ierr = PetscMalloc1(size+2,&owners_co);CHKERRQ(ierr);
1266f96d37fcSHong Zhang   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
1267f96d37fcSHong Zhang 
1268f96d37fcSHong Zhang   proc = 0;
1269f96d37fcSHong Zhang   for (i=0; i<pon; i++) {
1270f96d37fcSHong Zhang     while (prmap[i] >= owners[proc+1]) proc++;
1271f96d37fcSHong Zhang     len_si[proc]++;  /* num of rows in Co to be sent to [proc] */
1272f96d37fcSHong Zhang     len_s[proc] += coi[i+1] - coi[i];
1273f96d37fcSHong Zhang   }
1274f96d37fcSHong Zhang 
1275f96d37fcSHong Zhang   len          = 0; /* max length of buf_si[] */
1276f96d37fcSHong Zhang   owners_co[0] = 0;
1277f96d37fcSHong Zhang   for (proc=0; proc<size; proc++) {
1278f96d37fcSHong Zhang     owners_co[proc+1] = owners_co[proc] + len_si[proc];
1279f96d37fcSHong Zhang     if (len_si[proc]) {
1280f96d37fcSHong Zhang       merge->nsend++;
1281f96d37fcSHong Zhang       len_si[proc] = 2*(len_si[proc] + 1);
1282f96d37fcSHong Zhang       len         += len_si[proc];
1283f96d37fcSHong Zhang     }
1284f96d37fcSHong Zhang   }
1285f96d37fcSHong Zhang 
1286f96d37fcSHong Zhang   /* determine the number and length of messages to receive for coi and coj  */
12870298fd71SBarry Smith   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&merge->nrecv);CHKERRQ(ierr);
1288f96d37fcSHong Zhang   ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr);
1289f96d37fcSHong Zhang 
1290f96d37fcSHong Zhang   /* post the Irecv and Isend of coj */
1291f96d37fcSHong Zhang   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
1292f96d37fcSHong Zhang   ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
1293854ce69bSBarry Smith   ierr = PetscMalloc1(merge->nsend+1,&swaits);CHKERRQ(ierr);
1294f96d37fcSHong Zhang   for (proc=0, k=0; proc<size; proc++) {
1295f96d37fcSHong Zhang     if (!len_s[proc]) continue;
1296f96d37fcSHong Zhang     i    = owners_co[proc];
1297f96d37fcSHong Zhang     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
1298f96d37fcSHong Zhang     k++;
1299f96d37fcSHong Zhang   }
1300f96d37fcSHong Zhang 
1301f96d37fcSHong Zhang   /* receives and sends of coj are complete */
1302785e854fSJed Brown   ierr = PetscMalloc1(size,&sstatus);CHKERRQ(ierr);
1303f96d37fcSHong Zhang   for (i=0; i<merge->nrecv; i++) {
1304c280ed6aSJed Brown     PetscMPIInt icompleted;
1305c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1306f96d37fcSHong Zhang   }
1307f96d37fcSHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1308f96d37fcSHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1309f96d37fcSHong Zhang 
1310f96d37fcSHong Zhang   /* send and recv coi */
1311f96d37fcSHong Zhang   /*-------------------*/
1312f96d37fcSHong Zhang   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
1313f96d37fcSHong Zhang   ierr   = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
1314854ce69bSBarry Smith   ierr   = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr);
1315f96d37fcSHong Zhang   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
1316f96d37fcSHong Zhang   for (proc=0,k=0; proc<size; proc++) {
1317f96d37fcSHong Zhang     if (!len_s[proc]) continue;
1318f96d37fcSHong Zhang     /* form outgoing message for i-structure:
1319f96d37fcSHong Zhang          buf_si[0]:                 nrows to be sent
1320f96d37fcSHong Zhang                [1:nrows]:           row index (global)
1321f96d37fcSHong Zhang                [nrows+1:2*nrows+1]: i-structure index
1322f96d37fcSHong Zhang     */
1323f96d37fcSHong Zhang     /*-------------------------------------------*/
1324f96d37fcSHong Zhang     nrows       = len_si[proc]/2 - 1;
1325f96d37fcSHong Zhang     buf_si_i    = buf_si + nrows+1;
1326f96d37fcSHong Zhang     buf_si[0]   = nrows;
1327f96d37fcSHong Zhang     buf_si_i[0] = 0;
1328f96d37fcSHong Zhang     nrows       = 0;
1329f96d37fcSHong Zhang     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
1330f96d37fcSHong Zhang       nzi               = coi[i+1] - coi[i];
1331f96d37fcSHong Zhang       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */
1332f96d37fcSHong Zhang       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
1333f96d37fcSHong Zhang       nrows++;
1334f96d37fcSHong Zhang     }
1335f96d37fcSHong Zhang     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
1336f96d37fcSHong Zhang     k++;
1337f96d37fcSHong Zhang     buf_si += len_si[proc];
1338f96d37fcSHong Zhang   }
1339f96d37fcSHong Zhang   i = merge->nrecv;
1340f96d37fcSHong Zhang   while (i--) {
1341c280ed6aSJed Brown     PetscMPIInt icompleted;
1342c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1343f96d37fcSHong Zhang   }
1344f96d37fcSHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1345f96d37fcSHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1346f96d37fcSHong Zhang   ierr = PetscFree(len_si);CHKERRQ(ierr);
1347f96d37fcSHong Zhang   ierr = PetscFree(len_ri);CHKERRQ(ierr);
1348f96d37fcSHong Zhang   ierr = PetscFree(swaits);CHKERRQ(ierr);
1349f96d37fcSHong Zhang   ierr = PetscFree(sstatus);CHKERRQ(ierr);
1350f96d37fcSHong Zhang   ierr = PetscFree(buf_s);CHKERRQ(ierr);
1351f96d37fcSHong Zhang 
1352f96d37fcSHong Zhang   /* compute the local portion of C (mpi mat) */
1353f96d37fcSHong Zhang   /*------------------------------------------*/
1354f96d37fcSHong Zhang   /* allocate bi array and free space for accumulating nonzero column info */
1355854ce69bSBarry Smith   ierr  = PetscMalloc1(pn+1,&bi);CHKERRQ(ierr);
1356f96d37fcSHong Zhang   bi[0] = 0;
1357f96d37fcSHong Zhang 
1358c36aecf5SHong Zhang   /* set initial free space to be fill*(nnz(P) + nnz(A)) */
1359d6ab1dc0SHong Zhang   nnz           = fill*(pdti[pn] + poti[pon] + ai[am]);
1360a2ea699eSBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
1361f96d37fcSHong Zhang   current_space = free_space;
1362f96d37fcSHong Zhang 
1363dcca6d9dSJed Brown   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
1364f96d37fcSHong Zhang   for (k=0; k<merge->nrecv; k++) {
1365f96d37fcSHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1366f96d37fcSHong Zhang     nrows       = *buf_ri_k[k];
1367f96d37fcSHong Zhang     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
1368f96d37fcSHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
1369f96d37fcSHong Zhang   }
1370f96d37fcSHong Zhang 
13711c7d5954SHong Zhang   ierr = MatPreallocateInitialize(comm,pn,A->cmap->n,dnz,onz);CHKERRQ(ierr);
1372f96d37fcSHong Zhang   rmax = 0;
1373f96d37fcSHong Zhang   for (i=0; i<pn; i++) {
1374f96d37fcSHong Zhang     /* add pdt[i,:]*AP into lnk */
1375f96d37fcSHong Zhang     pnz = pdti[i+1] - pdti[i];
1376f96d37fcSHong Zhang     ptJ = pdtj + pdti[i];
1377f96d37fcSHong Zhang     for (j=0; j<pnz; j++) {
1378f96d37fcSHong Zhang       row  = ptJ[j];  /* row of AP == col of Pt */
1379d6ab1dc0SHong Zhang       anz  = ai[row+1] - ai[row];
1380d6ab1dc0SHong Zhang       Jptr = aj + ai[row];
1381f96d37fcSHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
1382d6ab1dc0SHong Zhang       ierr = PetscLLCondensedAddSorted(anz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1383f96d37fcSHong Zhang     }
13844e938277SHong Zhang 
1385f96d37fcSHong Zhang     /* add received col data into lnk */
1386f96d37fcSHong Zhang     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1387f96d37fcSHong Zhang       if (i == *nextrow[k]) { /* i-th row */
1388f96d37fcSHong Zhang         nzi  = *(nextci[k]+1) - *nextci[k];
1389f96d37fcSHong Zhang         Jptr = buf_rj[k] + *nextci[k];
13904e938277SHong Zhang         ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1391f96d37fcSHong Zhang         nextrow[k]++; nextci[k]++;
1392f96d37fcSHong Zhang       }
1393f96d37fcSHong Zhang     }
13944e938277SHong Zhang     nnz = lnk[0];
1395f96d37fcSHong Zhang 
1396f96d37fcSHong Zhang     /* if free space is not available, make more free space */
1397f96d37fcSHong Zhang     if (current_space->local_remaining<nnz) {
1398f96d37fcSHong Zhang       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
1399f96d37fcSHong Zhang       nspacedouble++;
1400f96d37fcSHong Zhang     }
1401f96d37fcSHong Zhang     /* copy data into free space, then initialize lnk */
14024e938277SHong Zhang     ierr = PetscLLCondensedClean(aN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
1403f96d37fcSHong Zhang     ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr);
14042205254eSKarl Rupp 
1405f96d37fcSHong Zhang     current_space->array           += nnz;
1406f96d37fcSHong Zhang     current_space->local_used      += nnz;
1407f96d37fcSHong Zhang     current_space->local_remaining -= nnz;
14082205254eSKarl Rupp 
1409f96d37fcSHong Zhang     bi[i+1] = bi[i] + nnz;
1410f96d37fcSHong Zhang     if (nnz > rmax) rmax = nnz;
1411f96d37fcSHong Zhang   }
1412f96d37fcSHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1413f96d37fcSHong Zhang 
1414854ce69bSBarry Smith   ierr = PetscMalloc1(bi[pn]+1,&bj);CHKERRQ(ierr);
1415f96d37fcSHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
14162205254eSKarl Rupp 
1417118727c9SMark F. Adams   afill_tmp = (PetscReal)bi[pn]/(pdti[pn] + poti[pon] + ai[am]+1);
1418f96d37fcSHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
1419d6ab1dc0SHong Zhang   ierr = PetscLLCondensedDestroy(lnk,lnkbt);CHKERRQ(ierr);
14204e938277SHong Zhang   ierr = MatDestroy(&POt);CHKERRQ(ierr);
14214e938277SHong Zhang   ierr = MatDestroy(&PDt);CHKERRQ(ierr);
1422f96d37fcSHong Zhang 
14231c7d5954SHong Zhang   /* create symbolic parallel matrix Cmpi - why cannot be assembled in Numeric part   */
14241c7d5954SHong Zhang   /*----------------------------------------------------------------------------------*/
14251795a4d1SJed Brown   ierr = PetscCalloc1(rmax+1,&vals);CHKERRQ(ierr);
1426f96d37fcSHong Zhang 
1427f96d37fcSHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
1428f96d37fcSHong Zhang   ierr = MatSetSizes(Cmpi,pn,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
142933d57670SJed Brown   ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(A->cmap->bs));CHKERRQ(ierr);
1430f96d37fcSHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
1431f96d37fcSHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
1432f96d37fcSHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
1433f96d37fcSHong Zhang   ierr = MatSetBlockSize(Cmpi,1);CHKERRQ(ierr);
1434f96d37fcSHong Zhang   for (i=0; i<pn; i++) {
1435f96d37fcSHong Zhang     row  = i + rstart;
1436f96d37fcSHong Zhang     nnz  = bi[i+1] - bi[i];
1437f96d37fcSHong Zhang     Jptr = bj + bi[i];
1438f96d37fcSHong Zhang     ierr = MatSetValues(Cmpi,1,&row,nnz,Jptr,vals,INSERT_VALUES);CHKERRQ(ierr);
1439f96d37fcSHong Zhang   }
1440f96d37fcSHong Zhang   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1441f96d37fcSHong Zhang   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1442f96d37fcSHong Zhang   ierr = PetscFree(vals);CHKERRQ(ierr);
1443f96d37fcSHong Zhang 
1444f96d37fcSHong Zhang   merge->bi        = bi;
1445f96d37fcSHong Zhang   merge->bj        = bj;
1446f96d37fcSHong Zhang   merge->coi       = coi;
1447f96d37fcSHong Zhang   merge->coj       = coj;
1448f96d37fcSHong Zhang   merge->buf_ri    = buf_ri;
1449f96d37fcSHong Zhang   merge->buf_rj    = buf_rj;
1450f96d37fcSHong Zhang   merge->owners_co = owners_co;
1451f96d37fcSHong Zhang   merge->destroy   = Cmpi->ops->destroy;
1452f96d37fcSHong Zhang   merge->duplicate = Cmpi->ops->duplicate;
1453f96d37fcSHong Zhang 
14546da69ca6SHong Zhang   Cmpi->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable;
1455c5af039cSHong Zhang   Cmpi->ops->destroy                 = MatDestroy_MPIAIJ_PtAP;
1456c9ba551fSBarry Smith   Cmpi->ops->duplicate               = MatDuplicate_MPIAIJ_MatPtAP;
1457f96d37fcSHong Zhang 
1458f96d37fcSHong Zhang   /* attach the supporting struct to Cmpi for reuse */
1459f96d37fcSHong Zhang   c           = (Mat_MPIAIJ*)Cmpi->data;
1460f96d37fcSHong Zhang   c->ptap     = ptap;
14610298fd71SBarry Smith   ptap->api   = NULL;
14620298fd71SBarry Smith   ptap->apj   = NULL;
1463f96d37fcSHong Zhang   ptap->merge = merge;
1464d6ab1dc0SHong Zhang   ptap->rmax  = rmax;
1465d6ab1dc0SHong Zhang 
1466d6ab1dc0SHong Zhang   *C = Cmpi;
1467d6ab1dc0SHong Zhang #if defined(PETSC_USE_INFO)
1468d6ab1dc0SHong Zhang   if (bi[pn] != 0) {
146957622a8eSBarry Smith     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
147057622a8eSBarry Smith     ierr = PetscInfo1(Cmpi,"Use MatTransposeMatMult(A,B,MatReuse,%g,&C) for best performance.\n",(double)afill);CHKERRQ(ierr);
1471d6ab1dc0SHong Zhang   } else {
1472d6ab1dc0SHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
1473d6ab1dc0SHong Zhang   }
1474d6ab1dc0SHong Zhang #endif
1475d6ab1dc0SHong Zhang   PetscFunctionReturn(0);
1476d6ab1dc0SHong Zhang }
1477d6ab1dc0SHong Zhang 
1478d6ab1dc0SHong Zhang #undef __FUNCT__
14796da69ca6SHong Zhang #define __FUNCT__ "MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ"
14806da69ca6SHong Zhang PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ(Mat P,Mat A,Mat C)
1481d6ab1dc0SHong Zhang {
1482d6ab1dc0SHong Zhang   PetscErrorCode      ierr;
1483d6ab1dc0SHong Zhang   Mat_Merge_SeqsToMPI *merge;
1484d6ab1dc0SHong Zhang   Mat_MPIAIJ          *p =(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
1485d6ab1dc0SHong Zhang   Mat_SeqAIJ          *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data;
1486d6ab1dc0SHong Zhang   Mat_PtAPMPI         *ptap;
1487e745005bSHong Zhang   PetscInt            *adj;
1488e745005bSHong Zhang   PetscInt            i,j,k,anz,pnz,row,*cj,nexta;
1489e745005bSHong Zhang   MatScalar           *ada,*ca,valtmp;
1490d6ab1dc0SHong Zhang   PetscInt            am  =A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n;
1491ce94432eSBarry Smith   MPI_Comm            comm;
1492d6ab1dc0SHong Zhang   PetscMPIInt         size,rank,taga,*len_s;
1493d6ab1dc0SHong Zhang   PetscInt            *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci;
1494d6ab1dc0SHong Zhang   PetscInt            **buf_ri,**buf_rj;
1495d6ab1dc0SHong Zhang   PetscInt            cnz=0,*bj_i,*bi,*bj,bnz,nextcj;  /* bi,bj,ba: local array of C(mpi mat) */
1496d6ab1dc0SHong Zhang   MPI_Request         *s_waits,*r_waits;
1497d6ab1dc0SHong Zhang   MPI_Status          *status;
1498d6ab1dc0SHong Zhang   MatScalar           **abuf_r,*ba_i,*pA,*coa,*ba;
1499d6ab1dc0SHong Zhang   PetscInt            *ai,*aj,*coi,*coj;
1500a2ea699eSBarry Smith   PetscInt            *poJ,*pdJ;
1501d6ab1dc0SHong Zhang   Mat                 A_loc;
1502d6ab1dc0SHong Zhang   Mat_SeqAIJ          *a_loc;
1503d6ab1dc0SHong Zhang 
1504d6ab1dc0SHong Zhang   PetscFunctionBegin;
1505ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr);
1506d6ab1dc0SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1507d6ab1dc0SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1508d6ab1dc0SHong Zhang 
1509d6ab1dc0SHong Zhang   ptap  = c->ptap;
1510d6ab1dc0SHong Zhang   merge = ptap->merge;
1511d6ab1dc0SHong Zhang 
1512e745005bSHong Zhang   /* 2) compute numeric C_seq = P_loc^T*A_loc */
1513e745005bSHong Zhang   /*------------------------------------------*/
1514d6ab1dc0SHong Zhang   /* get data from symbolic products */
1515d6ab1dc0SHong Zhang   coi    = merge->coi; coj = merge->coj;
1516854ce69bSBarry Smith   ierr   = PetscCalloc1(coi[pon]+1,&coa);CHKERRQ(ierr);
1517d6ab1dc0SHong Zhang   bi     = merge->bi; bj = merge->bj;
1518d6ab1dc0SHong Zhang   owners = merge->rowmap->range;
15191795a4d1SJed Brown   ierr   = PetscCalloc1(bi[cm]+1,&ba);CHKERRQ(ierr);
1520d6ab1dc0SHong Zhang 
1521d6ab1dc0SHong Zhang   /* get A_loc by taking all local rows of A */
1522d6ab1dc0SHong Zhang   A_loc = ptap->A_loc;
1523d6ab1dc0SHong Zhang   ierr  = MatMPIAIJGetLocalMat(A,MAT_REUSE_MATRIX,&A_loc);CHKERRQ(ierr);
1524d6ab1dc0SHong Zhang   a_loc = (Mat_SeqAIJ*)(A_loc)->data;
1525d6ab1dc0SHong Zhang   ai    = a_loc->i;
1526d6ab1dc0SHong Zhang   aj    = a_loc->j;
1527d6ab1dc0SHong Zhang 
1528d6ab1dc0SHong Zhang   for (i=0; i<am; i++) {
1529d6ab1dc0SHong Zhang     anz = ai[i+1] - ai[i];
1530d6ab1dc0SHong Zhang     adj = aj + ai[i];
1531d6ab1dc0SHong Zhang     ada = a_loc->a + ai[i];
1532d6ab1dc0SHong Zhang 
1533d6ab1dc0SHong Zhang     /* 2-b) Compute Cseq = P_loc[i,:]^T*A[i,:] using outer product */
1534e745005bSHong Zhang     /*-------------------------------------------------------------*/
1535d6ab1dc0SHong Zhang     /* put the value into Co=(p->B)^T*A (off-diagonal part, send to others) */
1536d6ab1dc0SHong Zhang     pnz = po->i[i+1] - po->i[i];
1537d6ab1dc0SHong Zhang     poJ = po->j + po->i[i];
1538d6ab1dc0SHong Zhang     pA  = po->a + po->i[i];
1539d6ab1dc0SHong Zhang     for (j=0; j<pnz; j++) {
1540d6ab1dc0SHong Zhang       row = poJ[j];
1541d6ab1dc0SHong Zhang       cj  = coj + coi[row];
1542d6ab1dc0SHong Zhang       ca  = coa + coi[row];
1543e745005bSHong Zhang       /* perform sparse axpy */
1544e745005bSHong Zhang       nexta  = 0;
1545d6ab1dc0SHong Zhang       valtmp = pA[j];
1546e745005bSHong Zhang       for (k=0; nexta<anz; k++) {
1547e745005bSHong Zhang         if (cj[k] == adj[nexta]) {
1548e745005bSHong Zhang           ca[k] += valtmp*ada[nexta];
1549e745005bSHong Zhang           nexta++;
1550d6ab1dc0SHong Zhang         }
1551e745005bSHong Zhang       }
1552e745005bSHong Zhang       ierr = PetscLogFlops(2.0*anz);CHKERRQ(ierr);
1553d6ab1dc0SHong Zhang     }
1554d6ab1dc0SHong Zhang 
1555d6ab1dc0SHong Zhang     /* put the value into Cd (diagonal part) */
1556d6ab1dc0SHong Zhang     pnz = pd->i[i+1] - pd->i[i];
1557d6ab1dc0SHong Zhang     pdJ = pd->j + pd->i[i];
1558d6ab1dc0SHong Zhang     pA  = pd->a + pd->i[i];
1559d6ab1dc0SHong Zhang     for (j=0; j<pnz; j++) {
1560d6ab1dc0SHong Zhang       row = pdJ[j];
1561d6ab1dc0SHong Zhang       cj  = bj + bi[row];
1562d6ab1dc0SHong Zhang       ca  = ba + bi[row];
1563e745005bSHong Zhang       /* perform sparse axpy */
1564e745005bSHong Zhang       nexta  = 0;
1565d6ab1dc0SHong Zhang       valtmp = pA[j];
1566e745005bSHong Zhang       for (k=0; nexta<anz; k++) {
1567e745005bSHong Zhang         if (cj[k] == adj[nexta]) {
1568e745005bSHong Zhang           ca[k] += valtmp*ada[nexta];
1569e745005bSHong Zhang           nexta++;
1570d6ab1dc0SHong Zhang         }
1571e745005bSHong Zhang       }
1572e745005bSHong Zhang       ierr = PetscLogFlops(2.0*anz);CHKERRQ(ierr);
1573d6ab1dc0SHong Zhang     }
1574d6ab1dc0SHong Zhang   }
1575d6ab1dc0SHong Zhang 
1576d6ab1dc0SHong Zhang   /* 3) send and recv matrix values coa */
1577d6ab1dc0SHong Zhang   /*------------------------------------*/
1578d6ab1dc0SHong Zhang   buf_ri = merge->buf_ri;
1579d6ab1dc0SHong Zhang   buf_rj = merge->buf_rj;
1580d6ab1dc0SHong Zhang   len_s  = merge->len_s;
1581d6ab1dc0SHong Zhang   ierr   = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr);
1582d6ab1dc0SHong Zhang   ierr   = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr);
1583d6ab1dc0SHong Zhang 
1584dcca6d9dSJed Brown   ierr = PetscMalloc2(merge->nsend+1,&s_waits,size,&status);CHKERRQ(ierr);
1585d6ab1dc0SHong Zhang   for (proc=0,k=0; proc<size; proc++) {
1586d6ab1dc0SHong Zhang     if (!len_s[proc]) continue;
1587d6ab1dc0SHong Zhang     i    = merge->owners_co[proc];
1588d6ab1dc0SHong Zhang     ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr);
1589d6ab1dc0SHong Zhang     k++;
1590d6ab1dc0SHong Zhang   }
1591d6ab1dc0SHong Zhang   if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);}
1592d6ab1dc0SHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);}
1593d6ab1dc0SHong Zhang 
1594d6ab1dc0SHong Zhang   ierr = PetscFree2(s_waits,status);CHKERRQ(ierr);
1595d6ab1dc0SHong Zhang   ierr = PetscFree(r_waits);CHKERRQ(ierr);
1596d6ab1dc0SHong Zhang   ierr = PetscFree(coa);CHKERRQ(ierr);
1597d6ab1dc0SHong Zhang 
1598d6ab1dc0SHong Zhang   /* 4) insert local Cseq and received values into Cmpi */
1599d6ab1dc0SHong Zhang   /*----------------------------------------------------*/
1600dcca6d9dSJed Brown   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
1601d6ab1dc0SHong Zhang   for (k=0; k<merge->nrecv; k++) {
1602e745005bSHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1603d6ab1dc0SHong Zhang     nrows       = *(buf_ri_k[k]);
1604d6ab1dc0SHong Zhang     nextrow[k]  = buf_ri_k[k]+1;  /* next row number of k-th recved i-structure */
1605d6ab1dc0SHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
1606d6ab1dc0SHong Zhang   }
1607d6ab1dc0SHong Zhang 
1608d6ab1dc0SHong Zhang   for (i=0; i<cm; i++) {
1609d6ab1dc0SHong Zhang     row  = owners[rank] + i; /* global row index of C_seq */
1610d6ab1dc0SHong Zhang     bj_i = bj + bi[i];  /* col indices of the i-th row of C */
1611d6ab1dc0SHong Zhang     ba_i = ba + bi[i];
1612d6ab1dc0SHong Zhang     bnz  = bi[i+1] - bi[i];
1613d6ab1dc0SHong Zhang     /* add received vals into ba */
1614d6ab1dc0SHong Zhang     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1615d6ab1dc0SHong Zhang       /* i-th row */
1616d6ab1dc0SHong Zhang       if (i == *nextrow[k]) {
1617d6ab1dc0SHong Zhang         cnz    = *(nextci[k]+1) - *nextci[k];
1618d6ab1dc0SHong Zhang         cj     = buf_rj[k] + *(nextci[k]);
1619d6ab1dc0SHong Zhang         ca     = abuf_r[k] + *(nextci[k]);
1620d6ab1dc0SHong Zhang         nextcj = 0;
1621d6ab1dc0SHong Zhang         for (j=0; nextcj<cnz; j++) {
1622d6ab1dc0SHong Zhang           if (bj_i[j] == cj[nextcj]) { /* bcol == ccol */
1623d6ab1dc0SHong Zhang             ba_i[j] += ca[nextcj++];
1624d6ab1dc0SHong Zhang           }
1625d6ab1dc0SHong Zhang         }
1626d6ab1dc0SHong Zhang         nextrow[k]++; nextci[k]++;
1627d6ab1dc0SHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1628d6ab1dc0SHong Zhang       }
1629d6ab1dc0SHong Zhang     }
1630d6ab1dc0SHong Zhang     ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr);
1631d6ab1dc0SHong Zhang   }
1632d6ab1dc0SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1633d6ab1dc0SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1634d6ab1dc0SHong Zhang 
1635d6ab1dc0SHong Zhang   ierr = PetscFree(ba);CHKERRQ(ierr);
1636d6ab1dc0SHong Zhang   ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr);
1637d6ab1dc0SHong Zhang   ierr = PetscFree(abuf_r);CHKERRQ(ierr);
1638d6ab1dc0SHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1639d6ab1dc0SHong Zhang   PetscFunctionReturn(0);
1640d6ab1dc0SHong Zhang }
1641d6ab1dc0SHong Zhang 
1642c9ba551fSBarry Smith PetscErrorCode MatDuplicate_MPIAIJ_MatPtAP(Mat, MatDuplicateOption,Mat*);
16432bbb1c24SHong Zhang /* This routine is modified from MatPtAPSymbolic_MPIAIJ_MPIAIJ();
16442bbb1c24SHong Zhang    differ from MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable in using LLCondensedCreate_Scalable() */
1645d6ab1dc0SHong Zhang #undef __FUNCT__
16466da69ca6SHong Zhang #define __FUNCT__ "MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ"
16476da69ca6SHong Zhang PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(Mat P,Mat A,PetscReal fill,Mat *C)
1648d6ab1dc0SHong Zhang {
1649d6ab1dc0SHong Zhang   PetscErrorCode      ierr;
1650d6ab1dc0SHong Zhang   Mat                 Cmpi,A_loc,POt,PDt;
1651d6ab1dc0SHong Zhang   Mat_PtAPMPI         *ptap;
16520298fd71SBarry Smith   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
1653d6ab1dc0SHong Zhang   Mat_MPIAIJ          *p        =(Mat_MPIAIJ*)P->data,*c;
1654d6ab1dc0SHong Zhang   PetscInt            *pdti,*pdtj,*poti,*potj,*ptJ;
1655d6ab1dc0SHong Zhang   PetscInt            nnz;
1656d6ab1dc0SHong Zhang   PetscInt            *lnk,*owners_co,*coi,*coj,i,k,pnz,row;
1657d6ab1dc0SHong Zhang   PetscInt            am  =A->rmap->n,pn=P->cmap->n;
1658ce94432eSBarry Smith   MPI_Comm            comm;
1659d6ab1dc0SHong Zhang   PetscMPIInt         size,rank,tagi,tagj,*len_si,*len_s,*len_ri;
1660d6ab1dc0SHong Zhang   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
1661d6ab1dc0SHong Zhang   PetscInt            len,proc,*dnz,*onz,*owners;
1662d6ab1dc0SHong Zhang   PetscInt            nzi,*bi,*bj;
1663d6ab1dc0SHong Zhang   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
1664d6ab1dc0SHong Zhang   MPI_Request         *swaits,*rwaits;
1665d6ab1dc0SHong Zhang   MPI_Status          *sstatus,rstatus;
1666d6ab1dc0SHong Zhang   Mat_Merge_SeqsToMPI *merge;
1667d6ab1dc0SHong Zhang   PetscInt            *ai,*aj,*Jptr,anz,*prmap=p->garray,pon,nspacedouble=0,j;
1668d6ab1dc0SHong Zhang   PetscReal           afill  =1.0,afill_tmp;
1669c36aecf5SHong Zhang   PetscInt            rstart = P->cmap->rstart,rmax,aN=A->cmap->N,Crmax;
1670d6ab1dc0SHong Zhang   PetscScalar         *vals;
1671d6ab1dc0SHong Zhang   Mat_SeqAIJ          *a_loc, *pdt,*pot;
1672d6ab1dc0SHong Zhang 
1673d6ab1dc0SHong Zhang   PetscFunctionBegin;
1674ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
1675d6ab1dc0SHong Zhang   /* check if matrix local sizes are compatible */
1676d6ab1dc0SHong Zhang   if (A->rmap->rstart != P->rmap->rstart || A->rmap->rend != P->rmap->rend) {
1677d6ab1dc0SHong Zhang     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);
1678d6ab1dc0SHong Zhang   }
1679d6ab1dc0SHong Zhang 
1680d6ab1dc0SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1681d6ab1dc0SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1682d6ab1dc0SHong Zhang 
1683d6ab1dc0SHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
1684b00a9115SJed Brown   ierr = PetscNew(&ptap);CHKERRQ(ierr);
1685d6ab1dc0SHong Zhang 
1686d6ab1dc0SHong Zhang   /* get A_loc by taking all local rows of A */
1687d6ab1dc0SHong Zhang   ierr = MatMPIAIJGetLocalMat(A,MAT_INITIAL_MATRIX,&A_loc);CHKERRQ(ierr);
16882205254eSKarl Rupp 
1689d6ab1dc0SHong Zhang   ptap->A_loc = A_loc;
1690d6ab1dc0SHong Zhang   a_loc       = (Mat_SeqAIJ*)(A_loc)->data;
1691d6ab1dc0SHong Zhang   ai          = a_loc->i;
1692d6ab1dc0SHong Zhang   aj          = a_loc->j;
1693d6ab1dc0SHong Zhang 
1694d6ab1dc0SHong Zhang   /* determine symbolic Co=(p->B)^T*A - send to others */
1695d6ab1dc0SHong Zhang   /*----------------------------------------------------*/
1696d6ab1dc0SHong Zhang   ierr = MatTransposeSymbolic_SeqAIJ(p->A,&PDt);CHKERRQ(ierr);
1697d6ab1dc0SHong Zhang   pdt  = (Mat_SeqAIJ*)PDt->data;
1698d6ab1dc0SHong Zhang   pdti = pdt->i; pdtj = pdt->j;
1699d6ab1dc0SHong Zhang 
1700d6ab1dc0SHong Zhang   ierr = MatTransposeSymbolic_SeqAIJ(p->B,&POt);CHKERRQ(ierr);
1701d6ab1dc0SHong Zhang   pot  = (Mat_SeqAIJ*)POt->data;
1702d6ab1dc0SHong Zhang   poti = pot->i; potj = pot->j;
1703d6ab1dc0SHong Zhang 
1704d6ab1dc0SHong Zhang   /* then, compute symbolic Co = (p->B)^T*A */
1705d6ab1dc0SHong Zhang   pon    = (p->B)->cmap->n; /* total num of rows to be sent to other processors
1706d6ab1dc0SHong Zhang                          >= (num of nonzero rows of C_seq) - pn */
1707854ce69bSBarry Smith   ierr   = PetscMalloc1(pon+1,&coi);CHKERRQ(ierr);
1708d6ab1dc0SHong Zhang   coi[0] = 0;
1709d6ab1dc0SHong Zhang 
1710d6ab1dc0SHong Zhang   /* set initial free space to be fill*(nnz(p->B) + nnz(A)) */
1711d6ab1dc0SHong Zhang   nnz           = fill*(poti[pon] + ai[am]);
1712a2ea699eSBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
1713d6ab1dc0SHong Zhang   current_space = free_space;
171419f0e57aSHong Zhang 
1715d6ab1dc0SHong Zhang   /* create and initialize a linked list */
1716d6ab1dc0SHong Zhang   i     = PetscMax(pdt->rmax,pot->rmax);
1717c36aecf5SHong Zhang   Crmax = i*a_loc->rmax*size; /* non-scalable! */
1718d6ab1dc0SHong Zhang   if (!Crmax || Crmax > aN) Crmax = aN;
1719d6ab1dc0SHong Zhang   ierr = PetscLLCondensedCreate_Scalable(Crmax,&lnk);CHKERRQ(ierr);
1720d6ab1dc0SHong Zhang 
1721d6ab1dc0SHong Zhang   for (i=0; i<pon; i++) {
1722d6ab1dc0SHong Zhang     pnz = poti[i+1] - poti[i];
1723d6ab1dc0SHong Zhang     ptJ = potj + poti[i];
1724d6ab1dc0SHong Zhang     for (j=0; j<pnz; j++) {
1725d6ab1dc0SHong Zhang       row  = ptJ[j]; /* row of A_loc == col of Pot */
1726d6ab1dc0SHong Zhang       anz  = ai[row+1] - ai[row];
1727d6ab1dc0SHong Zhang       Jptr = aj + ai[row];
1728d6ab1dc0SHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
1729d6ab1dc0SHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(anz,Jptr,lnk);CHKERRQ(ierr);
1730d6ab1dc0SHong Zhang     }
1731d6ab1dc0SHong Zhang     nnz = lnk[0];
1732d6ab1dc0SHong Zhang 
1733d6ab1dc0SHong Zhang     /* If free space is not available, double the total space in the list */
1734d6ab1dc0SHong Zhang     if (current_space->local_remaining<nnz) {
1735d6ab1dc0SHong Zhang       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
1736d6ab1dc0SHong Zhang       nspacedouble++;
1737d6ab1dc0SHong Zhang     }
1738d6ab1dc0SHong Zhang 
1739d6ab1dc0SHong Zhang     /* Copy data into free space, and zero out denserows */
1740d6ab1dc0SHong Zhang     ierr = PetscLLCondensedClean_Scalable(nnz,current_space->array,lnk);CHKERRQ(ierr);
17412205254eSKarl Rupp 
1742d6ab1dc0SHong Zhang     current_space->array           += nnz;
1743d6ab1dc0SHong Zhang     current_space->local_used      += nnz;
1744d6ab1dc0SHong Zhang     current_space->local_remaining -= nnz;
17452205254eSKarl Rupp 
1746d6ab1dc0SHong Zhang     coi[i+1] = coi[i] + nnz;
1747d6ab1dc0SHong Zhang   }
1748d6ab1dc0SHong Zhang 
1749854ce69bSBarry Smith   ierr = PetscMalloc1(coi[pon]+1,&coj);CHKERRQ(ierr);
1750d6ab1dc0SHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr);
17512205254eSKarl Rupp 
1752118727c9SMark F. Adams   afill_tmp = (PetscReal)coi[pon]/(poti[pon] + ai[am]+1);
1753d6ab1dc0SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
1754d6ab1dc0SHong Zhang 
1755d6ab1dc0SHong Zhang   /* send j-array (coj) of Co to other processors */
1756d6ab1dc0SHong Zhang   /*----------------------------------------------*/
1757d6ab1dc0SHong Zhang   /* determine row ownership */
1758b00a9115SJed Brown   ierr = PetscNew(&merge);CHKERRQ(ierr);
1759d6ab1dc0SHong Zhang   ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr);
17602205254eSKarl Rupp 
1761d6ab1dc0SHong Zhang   merge->rowmap->n  = pn;
1762d6ab1dc0SHong Zhang   merge->rowmap->bs = 1;
17632205254eSKarl Rupp 
1764d6ab1dc0SHong Zhang   ierr   = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr);
1765d6ab1dc0SHong Zhang   owners = merge->rowmap->range;
1766d6ab1dc0SHong Zhang 
1767d6ab1dc0SHong Zhang   /* determine the number of messages to send, their lengths */
17681795a4d1SJed Brown   ierr = PetscCalloc1(size,&len_si);CHKERRQ(ierr);
1769785e854fSJed Brown   ierr = PetscMalloc1(size,&merge->len_s);CHKERRQ(ierr);
17702205254eSKarl Rupp 
1771d6ab1dc0SHong Zhang   len_s        = merge->len_s;
1772d6ab1dc0SHong Zhang   merge->nsend = 0;
1773d6ab1dc0SHong Zhang 
1774854ce69bSBarry Smith   ierr = PetscMalloc1(size+2,&owners_co);CHKERRQ(ierr);
1775d6ab1dc0SHong Zhang   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
1776d6ab1dc0SHong Zhang 
1777d6ab1dc0SHong Zhang   proc = 0;
1778d6ab1dc0SHong Zhang   for (i=0; i<pon; i++) {
1779d6ab1dc0SHong Zhang     while (prmap[i] >= owners[proc+1]) proc++;
1780d6ab1dc0SHong Zhang     len_si[proc]++;  /* num of rows in Co to be sent to [proc] */
1781d6ab1dc0SHong Zhang     len_s[proc] += coi[i+1] - coi[i];
1782d6ab1dc0SHong Zhang   }
1783d6ab1dc0SHong Zhang 
1784d6ab1dc0SHong Zhang   len          = 0; /* max length of buf_si[] */
1785d6ab1dc0SHong Zhang   owners_co[0] = 0;
1786d6ab1dc0SHong Zhang   for (proc=0; proc<size; proc++) {
1787d6ab1dc0SHong Zhang     owners_co[proc+1] = owners_co[proc] + len_si[proc];
1788d6ab1dc0SHong Zhang     if (len_si[proc]) {
1789d6ab1dc0SHong Zhang       merge->nsend++;
1790d6ab1dc0SHong Zhang       len_si[proc] = 2*(len_si[proc] + 1);
1791d6ab1dc0SHong Zhang       len         += len_si[proc];
1792d6ab1dc0SHong Zhang     }
1793d6ab1dc0SHong Zhang   }
1794d6ab1dc0SHong Zhang 
1795d6ab1dc0SHong Zhang   /* determine the number and length of messages to receive for coi and coj  */
17960298fd71SBarry Smith   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&merge->nrecv);CHKERRQ(ierr);
1797d6ab1dc0SHong Zhang   ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr);
1798d6ab1dc0SHong Zhang 
1799d6ab1dc0SHong Zhang   /* post the Irecv and Isend of coj */
1800d6ab1dc0SHong Zhang   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
1801d6ab1dc0SHong Zhang   ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
1802854ce69bSBarry Smith   ierr = PetscMalloc1(merge->nsend+1,&swaits);CHKERRQ(ierr);
1803d6ab1dc0SHong Zhang   for (proc=0, k=0; proc<size; proc++) {
1804d6ab1dc0SHong Zhang     if (!len_s[proc]) continue;
1805d6ab1dc0SHong Zhang     i    = owners_co[proc];
1806d6ab1dc0SHong Zhang     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
1807d6ab1dc0SHong Zhang     k++;
1808d6ab1dc0SHong Zhang   }
1809d6ab1dc0SHong Zhang 
1810d6ab1dc0SHong Zhang   /* receives and sends of coj are complete */
1811785e854fSJed Brown   ierr = PetscMalloc1(size,&sstatus);CHKERRQ(ierr);
1812d6ab1dc0SHong Zhang   for (i=0; i<merge->nrecv; i++) {
1813c280ed6aSJed Brown     PetscMPIInt icompleted;
1814c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1815d6ab1dc0SHong Zhang   }
1816d6ab1dc0SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1817d6ab1dc0SHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1818d6ab1dc0SHong Zhang 
1819d6ab1dc0SHong Zhang   /* send and recv coi */
1820d6ab1dc0SHong Zhang   /*-------------------*/
1821d6ab1dc0SHong Zhang   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
1822d6ab1dc0SHong Zhang   ierr   = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
1823854ce69bSBarry Smith   ierr   = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr);
1824d6ab1dc0SHong Zhang   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
1825d6ab1dc0SHong Zhang   for (proc=0,k=0; proc<size; proc++) {
1826d6ab1dc0SHong Zhang     if (!len_s[proc]) continue;
1827d6ab1dc0SHong Zhang     /* form outgoing message for i-structure:
1828d6ab1dc0SHong Zhang          buf_si[0]:                 nrows to be sent
1829d6ab1dc0SHong Zhang                [1:nrows]:           row index (global)
1830d6ab1dc0SHong Zhang                [nrows+1:2*nrows+1]: i-structure index
1831d6ab1dc0SHong Zhang     */
1832d6ab1dc0SHong Zhang     /*-------------------------------------------*/
1833d6ab1dc0SHong Zhang     nrows       = len_si[proc]/2 - 1;
1834d6ab1dc0SHong Zhang     buf_si_i    = buf_si + nrows+1;
1835d6ab1dc0SHong Zhang     buf_si[0]   = nrows;
1836d6ab1dc0SHong Zhang     buf_si_i[0] = 0;
1837d6ab1dc0SHong Zhang     nrows       = 0;
1838d6ab1dc0SHong Zhang     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
1839d6ab1dc0SHong Zhang       nzi               = coi[i+1] - coi[i];
1840d6ab1dc0SHong Zhang       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi;  /* i-structure */
1841d6ab1dc0SHong Zhang       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
1842d6ab1dc0SHong Zhang       nrows++;
1843d6ab1dc0SHong Zhang     }
1844d6ab1dc0SHong Zhang     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
1845d6ab1dc0SHong Zhang     k++;
1846d6ab1dc0SHong Zhang     buf_si += len_si[proc];
1847d6ab1dc0SHong Zhang   }
1848d6ab1dc0SHong Zhang   i = merge->nrecv;
1849d6ab1dc0SHong Zhang   while (i--) {
1850c280ed6aSJed Brown     PetscMPIInt icompleted;
1851c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1852d6ab1dc0SHong Zhang   }
1853d6ab1dc0SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1854d6ab1dc0SHong Zhang   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1855d6ab1dc0SHong Zhang   ierr = PetscFree(len_si);CHKERRQ(ierr);
1856d6ab1dc0SHong Zhang   ierr = PetscFree(len_ri);CHKERRQ(ierr);
1857d6ab1dc0SHong Zhang   ierr = PetscFree(swaits);CHKERRQ(ierr);
1858d6ab1dc0SHong Zhang   ierr = PetscFree(sstatus);CHKERRQ(ierr);
1859d6ab1dc0SHong Zhang   ierr = PetscFree(buf_s);CHKERRQ(ierr);
1860d6ab1dc0SHong Zhang 
1861d6ab1dc0SHong Zhang   /* compute the local portion of C (mpi mat) */
1862d6ab1dc0SHong Zhang   /*------------------------------------------*/
1863d6ab1dc0SHong Zhang   /* allocate bi array and free space for accumulating nonzero column info */
1864854ce69bSBarry Smith   ierr  = PetscMalloc1(pn+1,&bi);CHKERRQ(ierr);
1865d6ab1dc0SHong Zhang   bi[0] = 0;
1866d6ab1dc0SHong Zhang 
1867d6ab1dc0SHong Zhang   /* set initial free space to be fill*(nnz(P) + nnz(AP)) */
1868d6ab1dc0SHong Zhang   nnz           = fill*(pdti[pn] + poti[pon] + ai[am]);
1869a2ea699eSBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
1870d6ab1dc0SHong Zhang   current_space = free_space;
1871d6ab1dc0SHong Zhang 
1872dcca6d9dSJed Brown   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
1873d6ab1dc0SHong Zhang   for (k=0; k<merge->nrecv; k++) {
1874d6ab1dc0SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1875d6ab1dc0SHong Zhang     nrows       = *buf_ri_k[k];
1876d6ab1dc0SHong Zhang     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
18772205254eSKarl Rupp     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* points to the next i-structure of k-th recieved i-structure  */
1878d6ab1dc0SHong Zhang   }
1879d6ab1dc0SHong Zhang 
1880d6ab1dc0SHong Zhang   ierr = MatPreallocateInitialize(comm,pn,A->cmap->n,dnz,onz);CHKERRQ(ierr);
1881d6ab1dc0SHong Zhang   rmax = 0;
1882d6ab1dc0SHong Zhang   for (i=0; i<pn; i++) {
1883d6ab1dc0SHong Zhang     /* add pdt[i,:]*AP into lnk */
1884d6ab1dc0SHong Zhang     pnz = pdti[i+1] - pdti[i];
1885d6ab1dc0SHong Zhang     ptJ = pdtj + pdti[i];
1886d6ab1dc0SHong Zhang     for (j=0; j<pnz; j++) {
1887d6ab1dc0SHong Zhang       row  = ptJ[j];  /* row of AP == col of Pt */
1888d6ab1dc0SHong Zhang       anz  = ai[row+1] - ai[row];
1889d6ab1dc0SHong Zhang       Jptr = aj + ai[row];
1890d6ab1dc0SHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
1891d6ab1dc0SHong Zhang       ierr = PetscLLCondensedAddSorted_Scalable(anz,Jptr,lnk);CHKERRQ(ierr);
1892d6ab1dc0SHong Zhang     }
1893d6ab1dc0SHong Zhang 
1894d6ab1dc0SHong Zhang     /* add received col data into lnk */
1895d6ab1dc0SHong Zhang     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1896d6ab1dc0SHong Zhang       if (i == *nextrow[k]) { /* i-th row */
1897d6ab1dc0SHong Zhang         nzi  = *(nextci[k]+1) - *nextci[k];
1898d6ab1dc0SHong Zhang         Jptr = buf_rj[k] + *nextci[k];
1899d6ab1dc0SHong Zhang         ierr = PetscLLCondensedAddSorted_Scalable(nzi,Jptr,lnk);CHKERRQ(ierr);
1900d6ab1dc0SHong Zhang         nextrow[k]++; nextci[k]++;
1901d6ab1dc0SHong Zhang       }
1902d6ab1dc0SHong Zhang     }
1903d6ab1dc0SHong Zhang     nnz = lnk[0];
1904d6ab1dc0SHong Zhang 
1905d6ab1dc0SHong Zhang     /* if free space is not available, make more free space */
1906d6ab1dc0SHong Zhang     if (current_space->local_remaining<nnz) {
1907d6ab1dc0SHong Zhang       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
1908d6ab1dc0SHong Zhang       nspacedouble++;
1909d6ab1dc0SHong Zhang     }
1910d6ab1dc0SHong Zhang     /* copy data into free space, then initialize lnk */
1911d6ab1dc0SHong Zhang     ierr = PetscLLCondensedClean_Scalable(nnz,current_space->array,lnk);CHKERRQ(ierr);
1912d6ab1dc0SHong Zhang     ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr);
19132205254eSKarl Rupp 
1914d6ab1dc0SHong Zhang     current_space->array           += nnz;
1915d6ab1dc0SHong Zhang     current_space->local_used      += nnz;
1916d6ab1dc0SHong Zhang     current_space->local_remaining -= nnz;
19172205254eSKarl Rupp 
1918d6ab1dc0SHong Zhang     bi[i+1] = bi[i] + nnz;
1919d6ab1dc0SHong Zhang     if (nnz > rmax) rmax = nnz;
1920d6ab1dc0SHong Zhang   }
19210d3441aeSHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr); //mpiexec -n 2 ./ex94 -f0 $D/tiny -f1 $D/tiny
1922d6ab1dc0SHong Zhang 
1923854ce69bSBarry Smith   ierr      = PetscMalloc1(bi[pn]+1,&bj);CHKERRQ(ierr);
1924d6ab1dc0SHong Zhang   ierr      = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
1925118727c9SMark F. Adams   afill_tmp = (PetscReal)bi[pn]/(pdti[pn] + poti[pon] + ai[am]+1);
1926d6ab1dc0SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
1927d6ab1dc0SHong Zhang   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr);
1928d6ab1dc0SHong Zhang   ierr = MatDestroy(&POt);CHKERRQ(ierr);
1929d6ab1dc0SHong Zhang   ierr = MatDestroy(&PDt);CHKERRQ(ierr);
1930d6ab1dc0SHong Zhang 
1931d6ab1dc0SHong Zhang   /* create symbolic parallel matrix Cmpi - why cannot be assembled in Numeric part   */
1932d6ab1dc0SHong Zhang   /*----------------------------------------------------------------------------------*/
19331795a4d1SJed Brown   ierr = PetscCalloc1(rmax+1,&vals);CHKERRQ(ierr);
1934d6ab1dc0SHong Zhang 
1935d6ab1dc0SHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
1936d6ab1dc0SHong Zhang   ierr = MatSetSizes(Cmpi,pn,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
193733d57670SJed Brown   ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(A->cmap->bs));CHKERRQ(ierr);
1938d6ab1dc0SHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
1939d6ab1dc0SHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
1940d6ab1dc0SHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
1941d6ab1dc0SHong Zhang   ierr = MatSetBlockSize(Cmpi,1);CHKERRQ(ierr);
1942d6ab1dc0SHong Zhang   for (i=0; i<pn; i++) {
1943d6ab1dc0SHong Zhang     row  = i + rstart;
1944d6ab1dc0SHong Zhang     nnz  = bi[i+1] - bi[i];
1945d6ab1dc0SHong Zhang     Jptr = bj + bi[i];
1946d6ab1dc0SHong Zhang     ierr = MatSetValues(Cmpi,1,&row,nnz,Jptr,vals,INSERT_VALUES);CHKERRQ(ierr);
1947d6ab1dc0SHong Zhang   }
1948d6ab1dc0SHong Zhang   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1949d6ab1dc0SHong Zhang   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1950d6ab1dc0SHong Zhang   ierr = PetscFree(vals);CHKERRQ(ierr);
1951d6ab1dc0SHong Zhang 
1952d6ab1dc0SHong Zhang   merge->bi        = bi;
1953d6ab1dc0SHong Zhang   merge->bj        = bj;
1954d6ab1dc0SHong Zhang   merge->coi       = coi;
1955d6ab1dc0SHong Zhang   merge->coj       = coj;
1956d6ab1dc0SHong Zhang   merge->buf_ri    = buf_ri;
1957d6ab1dc0SHong Zhang   merge->buf_rj    = buf_rj;
1958d6ab1dc0SHong Zhang   merge->owners_co = owners_co;
1959d6ab1dc0SHong Zhang   merge->destroy   = Cmpi->ops->destroy;
1960d6ab1dc0SHong Zhang   merge->duplicate = Cmpi->ops->duplicate;
1961d6ab1dc0SHong Zhang 
19626da69ca6SHong Zhang   Cmpi->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ;
1963d6ab1dc0SHong Zhang   Cmpi->ops->destroy                 = MatDestroy_MPIAIJ_PtAP;
1964c9ba551fSBarry Smith   Cmpi->ops->duplicate               = MatDuplicate_MPIAIJ_MatPtAP;
1965d6ab1dc0SHong Zhang 
1966d6ab1dc0SHong Zhang   /* attach the supporting struct to Cmpi for reuse */
1967d6ab1dc0SHong Zhang   c = (Mat_MPIAIJ*)Cmpi->data;
19682205254eSKarl Rupp 
1969d6ab1dc0SHong Zhang   c->ptap     = ptap;
19700298fd71SBarry Smith   ptap->api   = NULL;
19710298fd71SBarry Smith   ptap->apj   = NULL;
1972d6ab1dc0SHong Zhang   ptap->merge = merge;
1973d6ab1dc0SHong Zhang   ptap->rmax  = rmax;
19740298fd71SBarry Smith   ptap->apa   = NULL;
1975f96d37fcSHong Zhang 
1976f96d37fcSHong Zhang   *C = Cmpi;
1977f96d37fcSHong Zhang #if defined(PETSC_USE_INFO)
1978f96d37fcSHong Zhang   if (bi[pn] != 0) {
197957622a8eSBarry Smith     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
198057622a8eSBarry Smith     ierr = PetscInfo1(Cmpi,"Use MatTransposeMatMult(A,B,MatReuse,%g,&C) for best performance.\n",(double)afill);CHKERRQ(ierr);
1981f96d37fcSHong Zhang   } else {
1982f96d37fcSHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
1983f96d37fcSHong Zhang   }
1984f96d37fcSHong Zhang #endif
1985f96d37fcSHong Zhang   PetscFunctionReturn(0);
1986f96d37fcSHong Zhang }
1987