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