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