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