1 2 /* 3 Defines projective product routines where A is a MPIAIJ matrix 4 C = P^T * A * P 5 */ 6 7 #include <../src/mat/impls/aij/seq/aij.h> /*I "petscmat.h" I*/ 8 #include <../src/mat/utils/freespace.h> 9 #include <../src/mat/impls/aij/mpi/mpiaij.h> 10 #include <petscbt.h> 11 #include <petsctime.h> 12 13 #define PTAP_PROFILE 14 15 extern PetscErrorCode MatDestroy_MPIAIJ(Mat); 16 #undef __FUNCT__ 17 #define __FUNCT__ "MatDestroy_MPIAIJ_PtAP" 18 PetscErrorCode MatDestroy_MPIAIJ_PtAP(Mat A) 19 { 20 PetscErrorCode ierr; 21 Mat_MPIAIJ *a =(Mat_MPIAIJ*)A->data; 22 Mat_PtAPMPI *ptap=a->ptap; 23 24 PetscFunctionBegin; 25 if (ptap) { 26 Mat_Merge_SeqsToMPI *merge=ptap->merge; 27 ierr = PetscFree2(ptap->startsj_s,ptap->startsj_r);CHKERRQ(ierr); 28 ierr = PetscFree(ptap->bufa);CHKERRQ(ierr); 29 ierr = MatDestroy(&ptap->P_loc);CHKERRQ(ierr); 30 ierr = MatDestroy(&ptap->P_oth);CHKERRQ(ierr); 31 ierr = MatDestroy(&ptap->A_loc);CHKERRQ(ierr); /* used by MatTransposeMatMult() */ 32 ierr = MatDestroy(&ptap->R);CHKERRQ(ierr); 33 if (ptap->api) {ierr = PetscFree(ptap->api);CHKERRQ(ierr);} 34 if (ptap->apj) {ierr = PetscFree(ptap->apj);CHKERRQ(ierr);} 35 if (ptap->apa) {ierr = PetscFree(ptap->apa);CHKERRQ(ierr);} 36 if (merge) { 37 ierr = PetscFree(merge->id_r);CHKERRQ(ierr); 38 ierr = PetscFree(merge->len_s);CHKERRQ(ierr); 39 ierr = PetscFree(merge->len_r);CHKERRQ(ierr); 40 ierr = PetscFree(merge->bi);CHKERRQ(ierr); 41 ierr = PetscFree(merge->bj);CHKERRQ(ierr); 42 ierr = PetscFree(merge->buf_ri[0]);CHKERRQ(ierr); 43 ierr = PetscFree(merge->buf_ri);CHKERRQ(ierr); 44 ierr = PetscFree(merge->buf_rj[0]);CHKERRQ(ierr); 45 ierr = PetscFree(merge->buf_rj);CHKERRQ(ierr); 46 ierr = PetscFree(merge->coi);CHKERRQ(ierr); 47 ierr = PetscFree(merge->coj);CHKERRQ(ierr); 48 ierr = PetscFree(merge->owners_co);CHKERRQ(ierr); 49 ierr = PetscLayoutDestroy(&merge->rowmap);CHKERRQ(ierr); 50 ierr = merge->destroy(A);CHKERRQ(ierr); 51 ierr = PetscFree(ptap->merge);CHKERRQ(ierr); 52 } 53 ierr = PetscFree(ptap);CHKERRQ(ierr); 54 } 55 PetscFunctionReturn(0); 56 } 57 58 #undef __FUNCT__ 59 #define __FUNCT__ "MatDuplicate_MPIAIJ_MatPtAP" 60 PetscErrorCode MatDuplicate_MPIAIJ_MatPtAP(Mat A, MatDuplicateOption op, Mat *M) 61 { 62 PetscErrorCode ierr; 63 Mat_MPIAIJ *a = (Mat_MPIAIJ*)A->data; 64 Mat_PtAPMPI *ptap = a->ptap; 65 Mat_Merge_SeqsToMPI *merge = ptap->merge; 66 67 PetscFunctionBegin; 68 ierr = (*merge->duplicate)(A,op,M);CHKERRQ(ierr); 69 70 (*M)->ops->destroy = merge->destroy; 71 (*M)->ops->duplicate = merge->duplicate; 72 PetscFunctionReturn(0); 73 } 74 75 #undef __FUNCT__ 76 #define __FUNCT__ "MatPtAP_MPIAIJ_MPIAIJ" 77 PetscErrorCode MatPtAP_MPIAIJ_MPIAIJ(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C) 78 { 79 PetscErrorCode ierr; 80 81 PetscFunctionBegin; 82 if (scall == MAT_INITIAL_MATRIX) { 83 ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 84 ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ(A,P,fill,C);CHKERRQ(ierr); 85 ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 86 } 87 ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr); 88 ierr = MatPtAPNumeric_MPIAIJ_MPIAIJ(A,P,*C);CHKERRQ(ierr); 89 ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr); 90 PetscFunctionReturn(0); 91 } 92 93 #undef __FUNCT__ 94 #define __FUNCT__ "MatPtAPSymbolic_MPIAIJ_MPIAIJ" 95 PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C) 96 { 97 PetscErrorCode ierr; 98 Mat Cmpi; 99 Mat_PtAPMPI *ptap; 100 PetscFreeSpaceList free_space=NULL,current_space=NULL; 101 Mat_MPIAIJ *a =(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c; 102 Mat_SeqAIJ *ad =(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data; 103 Mat_SeqAIJ *p_loc,*p_oth; 104 PetscInt *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pdti,*pdtj,*poti,*potj,*ptJ; 105 PetscInt *adi=ad->i,*aj,*aoi=ao->i,nnz; 106 PetscInt *lnk,*owners_co,*coi,*coj,i,k,pnz,row; 107 PetscInt am=A->rmap->n,pN=P->cmap->N,pm=P->rmap->n,pn=P->cmap->n; 108 PetscBT lnkbt; 109 MPI_Comm comm; 110 PetscMPIInt size,rank,tagi,tagj,*len_si,*len_s,*len_ri,icompleted=0; 111 PetscInt **buf_rj,**buf_ri,**buf_ri_k; 112 PetscInt len,proc,*dnz,*onz,*owners; 113 PetscInt nzi,*pti,*ptj; 114 PetscInt nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci; 115 MPI_Request *swaits,*rwaits; 116 MPI_Status *sstatus,rstatus; 117 Mat_Merge_SeqsToMPI *merge; 118 PetscInt *api,*apj,*Jptr,apnz,*prmap=p->garray,pon,nspacedouble=0,j,ap_rmax=0; 119 PetscReal afill=1.0,afill_tmp; 120 PetscInt rmax; 121 122 PetscFunctionBegin; 123 ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 124 125 /* check if matrix local sizes are compatible */ 126 if (A->rmap->rstart != P->rmap->rstart || A->rmap->rend != P->rmap->rend) { 127 SETERRQ4(comm,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, Arow (%D, %D) != Prow (%D,%D)",A->rmap->rstart,A->rmap->rend,P->rmap->rstart,P->rmap->rend); 128 } 129 if (A->cmap->rstart != P->rmap->rstart || A->cmap->rend != P->rmap->rend) { 130 SETERRQ4(comm,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, Acol (%D, %D) != Prow (%D,%D)",A->cmap->rstart,A->cmap->rend,P->rmap->rstart,P->rmap->rend); 131 } 132 133 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 134 ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 135 136 /* create struct Mat_PtAPMPI and attached it to C later */ 137 ierr = PetscNew(&ptap);CHKERRQ(ierr); 138 ierr = PetscNew(&merge);CHKERRQ(ierr); 139 ptap->merge = merge; 140 ptap->reuse = MAT_INITIAL_MATRIX; 141 142 /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */ 143 ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 144 145 /* get P_loc by taking all local rows of P */ 146 ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 147 148 p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 149 p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 150 pi_loc = p_loc->i; pj_loc = p_loc->j; 151 pi_oth = p_oth->i; pj_oth = p_oth->j; 152 153 /* (1) compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth (api,apj) */ 154 /*--------------------------------------------------------------------------*/ 155 ierr = PetscMalloc1(am+1,&api);CHKERRQ(ierr); 156 api[0] = 0; 157 158 /* create and initialize a linked list */ 159 ierr = PetscLLCondensedCreate(pN,pN,&lnk,&lnkbt);CHKERRQ(ierr); 160 161 /* Initial FreeSpace size is fill*(nnz(A) + nnz(P)) -OOM for ex56, np=8k on Intrepid! */ 162 ierr = PetscFreeSpaceGet((PetscInt)(fill*(adi[am]+aoi[am]+pi_loc[pm])),&free_space);CHKERRQ(ierr); 163 164 current_space = free_space; 165 166 for (i=0; i<am; i++) { 167 /* diagonal portion of A */ 168 nzi = adi[i+1] - adi[i]; 169 aj = ad->j + adi[i]; 170 for (j=0; j<nzi; j++) { 171 row = aj[j]; 172 pnz = pi_loc[row+1] - pi_loc[row]; 173 Jptr = pj_loc + pi_loc[row]; 174 /* add non-zero cols of P into the sorted linked list lnk */ 175 ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 176 } 177 /* off-diagonal portion of A */ 178 nzi = aoi[i+1] - aoi[i]; 179 aj = ao->j + aoi[i]; 180 for (j=0; j<nzi; j++) { 181 row = aj[j]; 182 pnz = pi_oth[row+1] - pi_oth[row]; 183 Jptr = pj_oth + pi_oth[row]; 184 ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 185 } 186 apnz = lnk[0]; 187 api[i+1] = api[i] + apnz; 188 if (ap_rmax < apnz) ap_rmax = apnz; 189 190 /* if free space is not available, double the total space in the list */ 191 if (current_space->local_remaining<apnz) { 192 ierr = PetscFreeSpaceGet(apnz+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 193 nspacedouble++; 194 } 195 196 /* Copy data into free space, then initialize lnk */ 197 ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 198 199 current_space->array += apnz; 200 current_space->local_used += apnz; 201 current_space->local_remaining -= apnz; 202 } 203 204 /* Allocate space for apj, initialize apj, and */ 205 /* destroy list of free space and other temporary array(s) */ 206 ierr = PetscMalloc1(api[am]+1,&apj);CHKERRQ(ierr); 207 ierr = PetscFreeSpaceContiguous(&free_space,apj);CHKERRQ(ierr); 208 afill_tmp = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1); 209 if (afill_tmp > afill) afill = afill_tmp; 210 211 /* (2) determine symbolic Co=(p->B)^T*AP - send to others (coi,coj)*/ 212 /*-----------------------------------------------------------------*/ 213 ierr = MatGetSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr); 214 215 /* then, compute symbolic Co = (p->B)^T*AP */ 216 pon = (p->B)->cmap->n; /* total num of rows to be sent to other processors 217 >= (num of nonzero rows of C_seq) - pn */ 218 ierr = PetscMalloc1(pon+1,&coi);CHKERRQ(ierr); 219 coi[0] = 0; 220 221 /* set initial free space to be fill*(nnz(p->B) + nnz(AP)) */ 222 nnz = fill*(poti[pon] + api[am]); 223 ierr = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr); 224 current_space = free_space; 225 226 for (i=0; i<pon; i++) { 227 pnz = poti[i+1] - poti[i]; 228 ptJ = potj + poti[i]; 229 for (j=0; j<pnz; j++) { 230 row = ptJ[j]; /* row of AP == col of Pot */ 231 apnz = api[row+1] - api[row]; 232 Jptr = apj + api[row]; 233 /* add non-zero cols of AP into the sorted linked list lnk */ 234 ierr = PetscLLCondensedAddSorted(apnz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 235 } 236 nnz = lnk[0]; 237 238 /* If free space is not available, double the total space in the list */ 239 if (current_space->local_remaining<nnz) { 240 ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 241 nspacedouble++; 242 } 243 244 /* Copy data into free space, and zero out denserows */ 245 ierr = PetscLLCondensedClean(pN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 246 247 current_space->array += nnz; 248 current_space->local_used += nnz; 249 current_space->local_remaining -= nnz; 250 251 coi[i+1] = coi[i] + nnz; 252 } 253 254 ierr = PetscMalloc1(coi[pon],&coj);CHKERRQ(ierr); 255 ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr); 256 afill_tmp = (PetscReal)coi[pon]/(poti[pon] + api[am]+1); 257 if (afill_tmp > afill) afill = afill_tmp; 258 ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr); 259 260 /* (3) send j-array (coj) of Co to other processors */ 261 /*--------------------------------------------------*/ 262 ierr = PetscCalloc1(size,&merge->len_s);CHKERRQ(ierr); 263 len_s = merge->len_s; 264 merge->nsend = 0; 265 266 267 /* determine row ownership */ 268 ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr); 269 merge->rowmap->n = pn; 270 merge->rowmap->bs = 1; 271 272 ierr = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr); 273 owners = merge->rowmap->range; 274 275 /* determine the number of messages to send, their lengths */ 276 ierr = PetscMalloc2(size,&len_si,size,&sstatus);CHKERRQ(ierr); 277 ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 278 ierr = PetscMalloc1(size+2,&owners_co);CHKERRQ(ierr); 279 280 proc = 0; 281 for (i=0; i<pon; i++) { 282 while (prmap[i] >= owners[proc+1]) proc++; 283 len_si[proc]++; /* num of rows in Co(=Pt*AP) to be sent to [proc] */ 284 len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */ 285 } 286 287 len = 0; /* max length of buf_si[], see (4) */ 288 owners_co[0] = 0; 289 for (proc=0; proc<size; proc++) { 290 owners_co[proc+1] = owners_co[proc] + len_si[proc]; 291 if (len_s[proc]) { 292 merge->nsend++; 293 len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */ 294 len += len_si[proc]; 295 } 296 } 297 298 /* determine the number and length of messages to receive for coi and coj */ 299 ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&merge->nrecv);CHKERRQ(ierr); 300 ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr); 301 302 /* post the Irecv and Isend of coj */ 303 ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr); 304 ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr); 305 ierr = PetscMalloc1(merge->nsend+1,&swaits);CHKERRQ(ierr); 306 for (proc=0, k=0; proc<size; proc++) { 307 if (!len_s[proc]) continue; 308 i = owners_co[proc]; 309 ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr); 310 k++; 311 } 312 313 /* receives and sends of coj are complete */ 314 for (i=0; i<merge->nrecv; i++) { 315 ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 316 } 317 ierr = PetscFree(rwaits);CHKERRQ(ierr); 318 if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);} 319 320 /* (4) send and recv coi */ 321 /*-----------------------*/ 322 ierr = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr); 323 ierr = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr); 324 ierr = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr); 325 buf_si = buf_s; /* points to the beginning of k-th msg to be sent */ 326 for (proc=0,k=0; proc<size; proc++) { 327 if (!len_s[proc]) continue; 328 /* form outgoing message for i-structure: 329 buf_si[0]: nrows to be sent 330 [1:nrows]: row index (global) 331 [nrows+1:2*nrows+1]: i-structure index 332 */ 333 /*-------------------------------------------*/ 334 nrows = len_si[proc]/2 - 1; /* num of rows in Co to be sent to [proc] */ 335 buf_si_i = buf_si + nrows+1; 336 buf_si[0] = nrows; 337 buf_si_i[0] = 0; 338 nrows = 0; 339 for (i=owners_co[proc]; i<owners_co[proc+1]; i++) { 340 nzi = coi[i+1] - coi[i]; 341 buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */ 342 buf_si[nrows+1] = prmap[i] -owners[proc]; /* local row index */ 343 nrows++; 344 } 345 ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr); 346 k++; 347 buf_si += len_si[proc]; 348 } 349 i = merge->nrecv; 350 while (i--) { 351 ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 352 } 353 ierr = PetscFree(rwaits);CHKERRQ(ierr); 354 if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);} 355 356 ierr = PetscFree2(len_si,sstatus);CHKERRQ(ierr); 357 ierr = PetscFree(len_ri);CHKERRQ(ierr); 358 ierr = PetscFree(swaits);CHKERRQ(ierr); 359 ierr = PetscFree(buf_s);CHKERRQ(ierr); 360 361 /* (5) compute the local portion of C (mpi mat) */ 362 /*----------------------------------------------*/ 363 ierr = MatGetSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr); 364 365 /* allocate pti array and free space for accumulating nonzero column info */ 366 ierr = PetscMalloc1(pn+1,&pti);CHKERRQ(ierr); 367 pti[0] = 0; 368 369 /* set initial free space to be fill*(nnz(P) + nnz(AP)) */ 370 nnz = fill*(pi_loc[pm] + api[am]); 371 ierr = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr); 372 current_space = free_space; 373 374 ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr); 375 for (k=0; k<merge->nrecv; k++) { 376 buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 377 nrows = *buf_ri_k[k]; 378 nextrow[k] = buf_ri_k[k] + 1; /* next row number of k-th recved i-structure */ 379 nextci[k] = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure */ 380 } 381 ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr); 382 rmax = 0; 383 for (i=0; i<pn; i++) { 384 /* add pdt[i,:]*AP into lnk */ 385 pnz = pdti[i+1] - pdti[i]; 386 ptJ = pdtj + pdti[i]; 387 for (j=0; j<pnz; j++) { 388 row = ptJ[j]; /* row of AP == col of Pt */ 389 apnz = api[row+1] - api[row]; 390 Jptr = apj + api[row]; 391 /* add non-zero cols of AP into the sorted linked list lnk */ 392 ierr = PetscLLCondensedAddSorted(apnz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 393 } 394 395 /* add received col data into lnk */ 396 for (k=0; k<merge->nrecv; k++) { /* k-th received message */ 397 if (i == *nextrow[k]) { /* i-th row */ 398 nzi = *(nextci[k]+1) - *nextci[k]; 399 Jptr = buf_rj[k] + *nextci[k]; 400 ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr); 401 nextrow[k]++; nextci[k]++; 402 } 403 } 404 nnz = lnk[0]; 405 406 /* if free space is not available, make more free space */ 407 if (current_space->local_remaining<nnz) { 408 ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 409 nspacedouble++; 410 } 411 /* copy data into free space, then initialize lnk */ 412 ierr = PetscLLCondensedClean(pN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 413 ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr); 414 415 current_space->array += nnz; 416 current_space->local_used += nnz; 417 current_space->local_remaining -= nnz; 418 419 pti[i+1] = pti[i] + nnz; 420 if (nnz > rmax) rmax = nnz; 421 } 422 ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr); 423 ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr); 424 425 ierr = PetscMalloc1(pti[pn]+1,&ptj);CHKERRQ(ierr); 426 ierr = PetscFreeSpaceContiguous(&free_space,ptj);CHKERRQ(ierr); 427 afill_tmp = (PetscReal)pti[pn]/(pi_loc[pm] + api[am]+1); 428 if (afill_tmp > afill) afill = afill_tmp; 429 ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 430 431 /* (6) create symbolic parallel matrix Cmpi */ 432 /*------------------------------------------*/ 433 ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr); 434 ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 435 ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr); 436 ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr); 437 ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr); 438 ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 439 440 merge->bi = pti; /* Cseq->i */ 441 merge->bj = ptj; /* Cseq->j */ 442 merge->coi = coi; /* Co->i */ 443 merge->coj = coj; /* Co->j */ 444 merge->buf_ri = buf_ri; 445 merge->buf_rj = buf_rj; 446 merge->owners_co = owners_co; 447 merge->destroy = Cmpi->ops->destroy; 448 merge->duplicate = Cmpi->ops->duplicate; 449 450 /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */ 451 Cmpi->assembled = PETSC_FALSE; 452 Cmpi->ops->destroy = MatDestroy_MPIAIJ_PtAP; 453 Cmpi->ops->duplicate = MatDuplicate_MPIAIJ_MatPtAP; 454 455 /* attach the supporting struct to Cmpi for reuse */ 456 c = (Mat_MPIAIJ*)Cmpi->data; 457 c->ptap = ptap; 458 ptap->api = api; 459 ptap->apj = apj; 460 ptap->rmax = ap_rmax; 461 *C = Cmpi; 462 463 /* flag 'scalable' determines which implementations to be used: 464 0: do dense axpy in MatPtAPNumeric() - fast, but requires storage of a nonscalable dense array apa; 465 1: do sparse axpy in MatPtAPNumeric() - might slow, uses a sparse array apa */ 466 /* set default scalable */ 467 ptap->scalable = PETSC_FALSE; //PETSC_TRUE; 468 469 ierr = PetscOptionsGetBool(((PetscObject)Cmpi)->prefix,"-matptap_scalable",&ptap->scalable,NULL);CHKERRQ(ierr); 470 if (!ptap->scalable) { /* Do dense axpy */ 471 ierr = PetscCalloc1(pN,&ptap->apa);CHKERRQ(ierr); 472 ierr = MatTranspose_SeqAIJ(p->A,MAT_INITIAL_MATRIX,&ptap->R);CHKERRQ(ierr); 473 } else { 474 ierr = PetscCalloc1(ap_rmax+1,&ptap->apa);CHKERRQ(ierr); 475 } 476 477 #if defined(PETSC_USE_INFO) 478 if (pti[pn] != 0) { 479 ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr); 480 ierr = PetscInfo1(Cmpi,"Use MatPtAP(A,P,MatReuse,%g,&C) for best performance.\n",(double)afill);CHKERRQ(ierr); 481 } else { 482 ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr); 483 } 484 #endif 485 PetscFunctionReturn(0); 486 } 487 488 #undef __FUNCT__ 489 #define __FUNCT__ "MatPtAPNumeric_MPIAIJ_MPIAIJ" 490 PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C) 491 { 492 PetscErrorCode ierr; 493 Mat_MPIAIJ *a =(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data; 494 Mat_SeqAIJ *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data; 495 Mat_SeqAIJ *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data; 496 Mat_SeqAIJ *p_loc,*p_oth; 497 Mat_PtAPMPI *ptap; 498 Mat_Merge_SeqsToMPI *merge; 499 PetscInt *adi=ad->i,*aoi=ao->i,*adj,*aoj,*apJ,nextp; 500 PetscInt *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pJ,*pj; 501 PetscInt i,j,k,anz,pnz,apnz,nextap,row,*cj; 502 MatScalar *ada,*aoa,*apa,*pa,*ca,*pa_loc,*pa_oth,valtmp; 503 PetscInt am =A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n; 504 MPI_Comm comm; 505 PetscMPIInt size,rank,taga,*len_s; 506 PetscInt *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci; 507 PetscInt **buf_ri,**buf_rj; 508 PetscInt cnz=0,*bj_i,*bi,*bj,bnz,nextcj; /* bi,bj,ba: local array of C(mpi mat) */ 509 MPI_Request *s_waits,*r_waits; 510 MPI_Status *status; 511 MatScalar **abuf_r,*ba_i,*pA,*coa,*ba; 512 PetscInt *api,*apj,*coi,*coj; 513 PetscInt *poJ=po->j,*pdJ=pd->j,pcstart=P->cmap->rstart,pcend=P->cmap->rend; 514 PetscBool scalable; 515 #if defined(PTAP_PROFILE) 516 PetscLogDouble t0,t1,t2,t3,t4,et2_AP=0.0,et2_PtAP=0.0,t2_0,t2_1,t2_2,t_tran0,t_tran1; 517 #endif 518 519 PetscFunctionBegin; 520 ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr); 521 #if defined(PTAP_PROFILE) 522 ierr = PetscTime(&t0);CHKERRQ(ierr); 523 #endif 524 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 525 ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 526 527 ptap = c->ptap; 528 if (!ptap) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_ARG_INCOMP,"MatPtAP() has not been called to create matrix C yet, cannot use MAT_REUSE_MATRIX"); 529 merge = ptap->merge; 530 apa = ptap->apa; 531 scalable = ptap->scalable; 532 533 /* 1) get P_oth = ptap->P_oth and P_loc = ptap->P_loc */ 534 /*-----------------------------------------------------*/ 535 if (ptap->reuse == MAT_INITIAL_MATRIX) { 536 /* P_oth and P_loc are obtained in MatPtASymbolic(), skip calling MatGetBrowsOfAoCols() and MatMPIAIJGetLocalMat() */ 537 ptap->reuse = MAT_REUSE_MATRIX; 538 } else { /* update numerical values of P_oth and P_loc */ 539 ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 540 ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 541 } 542 #if defined(PTAP_PROFILE) 543 ierr = PetscTime(&t1);CHKERRQ(ierr); 544 #endif 545 546 /* 2) compute numeric C_seq = P_loc^T*A_loc*P - dominating part */ 547 /*--------------------------------------------------------------*/ 548 /* get data from symbolic products */ 549 p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 550 p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 551 pi_loc=p_loc->i; pj_loc=p_loc->j; pJ=pj_loc; pa_loc=p_loc->a; 552 pi_oth=p_oth->i; pj_oth=p_oth->j; pa_oth=p_oth->a; 553 554 coi = merge->coi; coj = merge->coj; 555 ierr = PetscCalloc1(coi[pon]+1,&coa);CHKERRQ(ierr); 556 557 bi = merge->bi; bj = merge->bj; 558 owners = merge->rowmap->range; 559 ierr = PetscCalloc1(bi[cm]+1,&ba);CHKERRQ(ierr); /* ba: Cseq->a */ 560 561 api = ptap->api; apj = ptap->apj; 562 563 if (!scalable) { /* Do dense axpy on apa (length of pN, stores A[i,:]*P) - nonscalable, but faster (could take 1/3 scalable time) */ 564 ierr = PetscInfo(C,"Using non-scalable dense axpy\n");CHKERRQ(ierr); 565 /*-----------------------------------------------------------------------------------------------------*/ 566 ierr = PetscTime(&t_tran0);CHKERRQ(ierr); 567 ierr = MatTranspose_SeqAIJ(p->A,MAT_REUSE_MATRIX,&ptap->R);CHKERRQ(ierr); 568 ierr = PetscTime(&t_tran1);CHKERRQ(ierr); 569 570 Mat R = ptap->R; 571 Mat_SeqAIJ *r = (Mat_SeqAIJ*)R->data; 572 PetscInt *ri=r->i,*rj=r->j,rnz,arow,l,prow,pcol,pN=P->cmap->N; 573 PetscScalar *ra=r->a,tmp,cdense[pN]; 574 575 ierr = PetscMemzero(cdense,pN*sizeof(PetscScalar));CHKERRQ(ierr); 576 for (i=0; i<cm; i++) { /* each row of C or R */ 577 rnz = ri[i+1] - ri[i]; 578 579 for (j=0; j<rnz; j++) { /* each nz of R */ 580 arow = rj[ri[i] + j]; 581 582 /* diagonal portion of A */ 583 anz = ad->i[arow+1] - ad->i[arow]; 584 for (k=0; k<anz; k++) { /* each nz of Ad */ 585 tmp = ra[ri[i] + j]*ad->a[ad->i[arow] + k]; 586 prow = ad->j[ad->i[arow] + k]; 587 pnz = pi_loc[prow+1] - pi_loc[prow]; 588 589 for (l=0; l<pnz; l++) { /* each nz of P_loc */ 590 pcol = pj_loc[pi_loc[prow] + l]; 591 cdense[pcol] += tmp*pa_loc[pi_loc[prow] + l]; 592 } 593 } 594 595 /* off-diagonal portion of A */ 596 anz = ao->i[arow+1] - ao->i[arow]; 597 for (k=0; k<anz; k++) { /* each nz of Ao */ 598 tmp = ra[ri[i] + j]*ao->a[ao->i[arow] + k]; 599 prow = ao->j[ao->i[arow] + k]; 600 pnz = pi_oth[prow+1] - pi_oth[prow]; 601 602 for (l=0; l<pnz; l++) { /* each nz of P_oth */ 603 pcol = pj_oth[pi_oth[prow] + l]; 604 cdense[pcol] += tmp*pa_oth[pi_oth[prow] + l]; 605 } 606 } 607 608 } //for (j=0; j<rnz; j++) 609 610 /* copy cdense[] into ca; zero cdense[] */ 611 cnz = bi[i+1] - bi[i]; 612 cj = bj + bi[i]; 613 ca = ba + bi[i]; 614 for (j=0; j<cnz; j++) { 615 ca[j] += cdense[cj[j]]; 616 cdense[cj[j]] = 0.0; 617 } 618 #if 0 619 if (rank == 0) { 620 printf("[%d] row %d: ",rank,i); 621 for (j=0; j<pN; j++) printf(" %g,",cdense[j]); 622 printf("\n"); 623 } 624 for (j=0; j<pN; j++) cdense[j]=0.0; // zero cdnese[] 625 #endif 626 } //for (i=0; i<cm; i++) { 627 628 //========================================== 629 630 631 for (i=0; i<am; i++) { 632 #if defined(PTAP_PROFILE) 633 ierr = PetscTime(&t2_0);CHKERRQ(ierr); 634 #endif 635 /* 2-a) form i-th sparse row of A_loc*P = Ad*P_loc + Ao*P_oth */ 636 /*------------------------------------------------------------*/ 637 apJ = apj + api[i]; 638 639 /* diagonal portion of A */ 640 anz = adi[i+1] - adi[i]; 641 adj = ad->j + adi[i]; 642 ada = ad->a + adi[i]; 643 for (j=0; j<anz; j++) { 644 row = adj[j]; 645 pnz = pi_loc[row+1] - pi_loc[row]; 646 pj = pj_loc + pi_loc[row]; 647 pa = pa_loc + pi_loc[row]; 648 649 /* perform dense axpy */ 650 valtmp = ada[j]; 651 for (k=0; k<pnz; k++) { 652 apa[pj[k]] += valtmp*pa[k]; 653 } 654 ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr); 655 } 656 657 /* off-diagonal portion of A */ 658 anz = aoi[i+1] - aoi[i]; 659 aoj = ao->j + aoi[i]; 660 aoa = ao->a + aoi[i]; 661 for (j=0; j<anz; j++) { 662 row = aoj[j]; 663 pnz = pi_oth[row+1] - pi_oth[row]; 664 pj = pj_oth + pi_oth[row]; 665 pa = pa_oth + pi_oth[row]; 666 667 /* perform dense axpy */ 668 valtmp = aoa[j]; 669 for (k=0; k<pnz; k++) { 670 apa[pj[k]] += valtmp*pa[k]; 671 } 672 ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr); 673 } 674 #if defined(PTAP_PROFILE) 675 ierr = PetscTime(&t2_1);CHKERRQ(ierr); 676 et2_AP += t2_1 - t2_0; 677 #endif 678 679 /* 2-b) Compute Cseq = P_loc[i,:]^T*AP[i,:] using outer product */ 680 /*--------------------------------------------------------------*/ 681 apnz = api[i+1] - api[i]; 682 /* put the value into Co=(p->B)^T*AP (off-diagonal part, send to others) */ 683 pnz = po->i[i+1] - po->i[i]; 684 poJ = po->j + po->i[i]; 685 pA = po->a + po->i[i]; 686 for (j=0; j<pnz; j++) { 687 row = poJ[j]; 688 cnz = coi[row+1] - coi[row]; 689 cj = coj + coi[row]; 690 ca = coa + coi[row]; 691 /* perform dense axpy */ 692 valtmp = pA[j]; 693 for (k=0; k<cnz; k++) { 694 ca[k] += valtmp*apa[cj[k]]; 695 } 696 ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr); 697 } 698 #if 0 699 /* put the value into Cd (diagonal part) */ 700 pnz = pd->i[i+1] - pd->i[i]; 701 pdJ = pd->j + pd->i[i]; 702 pA = pd->a + pd->i[i]; 703 for (j=0; j<pnz; j++) { 704 row = pdJ[j]; 705 cnz = bi[row+1] - bi[row]; 706 cj = bj + bi[row]; 707 ca = ba + bi[row]; 708 /* perform dense axpy */ 709 valtmp = pA[j]; 710 for (k=0; k<cnz; k++) { 711 ca[k] += valtmp*apa[cj[k]]; 712 } 713 ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr); 714 } 715 #endif 716 /* zero the current row of A*P */ 717 for (k=0; k<apnz; k++) apa[apJ[k]] = 0.0; 718 #if defined(PTAP_PROFILE) 719 ierr = PetscTime(&t2_2);CHKERRQ(ierr); 720 et2_PtAP += t2_2 - t2_1; 721 #endif 722 } 723 724 if (rank == 100) { 725 for (row=0; row<cm; row++) { 726 printf("[%d] row %d: ",rank,row); 727 cnz = bi[row+1] - bi[row]; 728 for (j=0; j<cnz; j++) printf(" %g,",ba[bi[row]+j]); 729 printf("\n"); 730 } 731 } 732 733 } else { /* Do sparse axpy on apa (length of ap_rmax, stores A[i,:]*P) - scalable, but slower */ 734 ierr = PetscInfo(C,"Using scalable sparse axpy\n");CHKERRQ(ierr); 735 /*-----------------------------------------------------------------------------------------*/ 736 pA=pa_loc; 737 for (i=0; i<am; i++) { 738 #if defined(PTAP_PROFILE) 739 ierr = PetscTime(&t2_0);CHKERRQ(ierr); 740 #endif 741 /* form i-th sparse row of A*P */ 742 apnz = api[i+1] - api[i]; 743 apJ = apj + api[i]; 744 /* diagonal portion of A */ 745 anz = adi[i+1] - adi[i]; 746 adj = ad->j + adi[i]; 747 ada = ad->a + adi[i]; 748 for (j=0; j<anz; j++) { 749 row = adj[j]; 750 pnz = pi_loc[row+1] - pi_loc[row]; 751 pj = pj_loc + pi_loc[row]; 752 pa = pa_loc + pi_loc[row]; 753 valtmp = ada[j]; 754 nextp = 0; 755 for (k=0; nextp<pnz; k++) { 756 if (apJ[k] == pj[nextp]) { /* col of AP == col of P */ 757 apa[k] += valtmp*pa[nextp++]; 758 } 759 } 760 ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr); 761 } 762 /* off-diagonal portion of A */ 763 anz = aoi[i+1] - aoi[i]; 764 aoj = ao->j + aoi[i]; 765 aoa = ao->a + aoi[i]; 766 for (j=0; j<anz; j++) { 767 row = aoj[j]; 768 pnz = pi_oth[row+1] - pi_oth[row]; 769 pj = pj_oth + pi_oth[row]; 770 pa = pa_oth + pi_oth[row]; 771 valtmp = aoa[j]; 772 nextp = 0; 773 for (k=0; nextp<pnz; k++) { 774 if (apJ[k] == pj[nextp]) { /* col of AP == col of P */ 775 apa[k] += valtmp*pa[nextp++]; 776 } 777 } 778 ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr); 779 } 780 #if defined(PTAP_PROFILE) 781 ierr = PetscTime(&t2_1);CHKERRQ(ierr); 782 et2_AP += t2_1 - t2_0; 783 #endif 784 785 /* 2-b) Compute Cseq = P_loc[i,:]^T*AP[i,:] using outer product */ 786 /*--------------------------------------------------------------*/ 787 pnz = pi_loc[i+1] - pi_loc[i]; 788 pJ = pj_loc + pi_loc[i]; 789 for (j=0; j<pnz; j++) { 790 nextap = 0; 791 row = pJ[j]; /* global index */ 792 if (row < pcstart || row >=pcend) { /* put the value into Co */ 793 row = *poJ; 794 cj = coj + coi[row]; 795 ca = coa + coi[row]; poJ++; 796 } else { /* put the value into Cd */ 797 row = *pdJ; 798 cj = bj + bi[row]; 799 ca = ba + bi[row]; pdJ++; 800 } 801 valtmp = pA[j]; 802 for (k=0; nextap<apnz; k++) { 803 if (cj[k]==apJ[nextap]) ca[k] += valtmp*apa[nextap++]; 804 } 805 ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr); 806 } 807 pA += pnz; 808 /* zero the current row info for A*P */ 809 ierr = PetscMemzero(apa,apnz*sizeof(MatScalar));CHKERRQ(ierr); 810 #if defined(PTAP_PROFILE) 811 ierr = PetscTime(&t2_2);CHKERRQ(ierr); 812 et2_PtAP += t2_2 - t2_1; 813 #endif 814 } 815 } 816 #if defined(PTAP_PROFILE) 817 ierr = PetscTime(&t2);CHKERRQ(ierr); 818 #endif 819 820 /* 3) send and recv matrix values coa */ 821 /*------------------------------------*/ 822 buf_ri = merge->buf_ri; 823 buf_rj = merge->buf_rj; 824 len_s = merge->len_s; 825 ierr = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr); 826 ierr = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr); 827 828 ierr = PetscMalloc2(merge->nsend+1,&s_waits,size,&status);CHKERRQ(ierr); 829 for (proc=0,k=0; proc<size; proc++) { 830 if (!len_s[proc]) continue; 831 i = merge->owners_co[proc]; 832 ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr); 833 k++; 834 } 835 if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);} 836 if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);} 837 838 ierr = PetscFree2(s_waits,status);CHKERRQ(ierr); 839 ierr = PetscFree(r_waits);CHKERRQ(ierr); 840 ierr = PetscFree(coa);CHKERRQ(ierr); 841 #if defined(PTAP_PROFILE) 842 ierr = PetscTime(&t3);CHKERRQ(ierr); 843 #endif 844 845 /* 4) insert local Cseq and received values into Cmpi */ 846 /*------------------------------------------------------*/ 847 ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr); 848 for (k=0; k<merge->nrecv; k++) { 849 buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 850 nrows = *(buf_ri_k[k]); 851 nextrow[k] = buf_ri_k[k]+1; /* next row number of k-th recved i-structure */ 852 nextci[k] = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure */ 853 } 854 855 for (i=0; i<cm; i++) { 856 row = owners[rank] + i; /* global row index of C_seq */ 857 bj_i = bj + bi[i]; /* col indices of the i-th row of C */ 858 ba_i = ba + bi[i]; 859 bnz = bi[i+1] - bi[i]; 860 /* add received vals into ba */ 861 for (k=0; k<merge->nrecv; k++) { /* k-th received message */ 862 /* i-th row */ 863 if (i == *nextrow[k]) { 864 cnz = *(nextci[k]+1) - *nextci[k]; 865 cj = buf_rj[k] + *(nextci[k]); 866 ca = abuf_r[k] + *(nextci[k]); 867 nextcj = 0; 868 for (j=0; nextcj<cnz; j++) { 869 if (bj_i[j] == cj[nextcj]) { /* bcol == ccol */ 870 ba_i[j] += ca[nextcj++]; 871 } 872 } 873 nextrow[k]++; nextci[k]++; 874 ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr); 875 } 876 } 877 ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr); 878 } 879 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 880 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 881 882 ierr = PetscFree(ba);CHKERRQ(ierr); 883 ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr); 884 ierr = PetscFree(abuf_r);CHKERRQ(ierr); 885 ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr); 886 #if defined(PTAP_PROFILE) 887 ierr = PetscTime(&t4);CHKERRQ(ierr); 888 if (rank==1) { 889 ierr = PetscPrintf(MPI_COMM_SELF," R=Pd^T %g\n", t_tran1 - t_tran0);CHKERRQ(ierr); 890 ierr = PetscPrintf(MPI_COMM_SELF," [%d] PtAPNum %g/P + %g/PtAP( %g/A*P + %g/Pt*AP ) + %g/comm + %g/Cloc = %g\n\n",rank,t1-t0,t2-t1,et2_AP,et2_PtAP,t3-t2,t4-t3,t4-t0);CHKERRQ(ierr); 891 } 892 #endif 893 PetscFunctionReturn(0); 894 } 895