1 #define PETSCMAT_DLL 2 3 /* 4 Defines projective product routines where A is a MPIAIJ matrix 5 C = P^T * A * P 6 */ 7 8 #include "src/mat/impls/aij/seq/aij.h" /*I "petscmat.h" I*/ 9 #include "src/mat/utils/freespace.h" 10 #include "src/mat/impls/aij/mpi/mpiaij.h" 11 #include "petscbt.h" 12 13 #undef __FUNCT__ 14 #define __FUNCT__ "MatPtAPSymbolic_MPIAIJ" 15 PetscErrorCode MatPtAPSymbolic_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C) 16 { 17 PetscErrorCode ierr; 18 19 PetscFunctionBegin; 20 if (!P->ops->ptapsymbolic_mpiaij) { 21 SETERRQ2(PETSC_ERR_SUP,"Not implemented for A=%s and P=%s",A->type_name,P->type_name); 22 } 23 ierr = (*P->ops->ptapsymbolic_mpiaij)(A,P,fill,C);CHKERRQ(ierr); 24 PetscFunctionReturn(0); 25 } 26 27 #undef __FUNCT__ 28 #define __FUNCT__ "MatPtAPNumeric_MPIAIJ" 29 PetscErrorCode MatPtAPNumeric_MPIAIJ(Mat A,Mat P,Mat C) 30 { 31 PetscErrorCode ierr; 32 33 PetscFunctionBegin; 34 if (!P->ops->ptapnumeric_mpiaij) { 35 SETERRQ2(PETSC_ERR_SUP,"Not implemented for A=%s and P=%s",A->type_name,P->type_name); 36 } 37 ierr = (*P->ops->ptapnumeric_mpiaij)(A,P,C);CHKERRQ(ierr); 38 PetscFunctionReturn(0); 39 } 40 41 #undef __FUNCT__ 42 #define __FUNCT__ "MatPtAPSymbolic_MPIAIJ_MPIAIJ" 43 PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C) 44 { 45 PetscErrorCode ierr; 46 Mat B_mpi; 47 Mat_MatMatMultMPI *ap; 48 PetscObjectContainer container; 49 FreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 50 Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data; 51 Mat_SeqAIJ *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data; 52 Mat_SeqAIJ *p_loc,*p_oth; 53 PetscInt *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pdti,*pdtj,*poti,*potj,*ptJ; 54 PetscInt *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,nnz; 55 PetscInt nlnk,*lnk,*owners_co,*coi,*coj,i,k,pnz,row; 56 PetscInt am=A->m,pN=P->N,pn=P->n; 57 PetscBT lnkbt; 58 MPI_Comm comm=A->comm; 59 PetscMPIInt size,rank,tag,*len_si,*len_s,*len_ri; 60 PetscInt **buf_rj,**buf_ri,**buf_ri_k; 61 PetscInt len,proc,*dnz,*onz,*owners; 62 PetscInt nzi,*bi,*bj; 63 PetscInt nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci; 64 MPI_Request *swaits,*rwaits; 65 MPI_Status *sstatus,rstatus; 66 Mat_Merge_SeqsToMPI *merge; 67 PetscInt *api,*apj,*Jptr,apnz,*prmap=p->garray,pon; 68 PetscMPIInt j; 69 70 PetscFunctionBegin; 71 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 72 ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 73 74 /* destroy the container 'Mat_MatMatMultMPI' in case that P is attached to it */ 75 ierr = PetscObjectQuery((PetscObject)P,"Mat_MatMatMultMPI",(PetscObject *)&container);CHKERRQ(ierr); 76 if (container) { 77 ierr = PetscObjectContainerDestroy(container);CHKERRQ(ierr); 78 ierr = PetscObjectCompose((PetscObject)P,"Mat_MatMatMultMPI",0);CHKERRQ(ierr); 79 } 80 81 /* create the container 'Mat_MatMatMultMPI' and attach it to P */ 82 ierr = PetscNew(Mat_MatMatMultMPI,&ap);CHKERRQ(ierr); 83 ap->abi=PETSC_NULL; ap->abj=PETSC_NULL; 84 ap->abnz_max = 0; 85 86 ierr = PetscObjectContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr); 87 ierr = PetscObjectContainerSetPointer(container,ap);CHKERRQ(ierr); 88 ierr = PetscObjectCompose((PetscObject)P,"Mat_MatMatMultMPI",(PetscObject)container);CHKERRQ(ierr); 89 ap->MatDestroy = P->ops->destroy; 90 P->ops->destroy = MatDestroy_MPIAIJ_MatMatMult; 91 ap->reuse = MAT_INITIAL_MATRIX; 92 ierr = PetscObjectContainerSetUserDestroy(container,PetscObjectContainerDestroy_Mat_MatMatMultMPI);CHKERRQ(ierr); 93 94 /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */ 95 ierr = MatGetBrowsOfAoCols(A,P,MAT_INITIAL_MATRIX,&ap->startsj,&ap->bufa,&ap->B_oth);CHKERRQ(ierr); 96 /* get P_loc by taking all local rows of P */ 97 ierr = MatGetLocalMat(P,MAT_INITIAL_MATRIX,&ap->B_loc);CHKERRQ(ierr); 98 99 p_loc = (Mat_SeqAIJ*)(ap->B_loc)->data; 100 p_oth = (Mat_SeqAIJ*)(ap->B_oth)->data; 101 pi_loc = p_loc->i; pj_loc = p_loc->j; 102 pi_oth = p_oth->i; pj_oth = p_oth->j; 103 104 /* first, compute symbolic AP = A_loc*P */ 105 /*--------------------------------------*/ 106 ierr = PetscMalloc((am+2)*sizeof(PetscInt),&api);CHKERRQ(ierr); 107 ap->abi = api; 108 api[0] = 0; 109 110 /* create and initialize a linked list */ 111 nlnk = pN+1; 112 ierr = PetscLLCreate(pN,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr); 113 114 /* Initial FreeSpace size is fill*nnz(A) */ 115 ierr = GetMoreSpace((PetscInt)(fill*(adi[am]+aoi[am])),&free_space);CHKERRQ(ierr); 116 current_space = free_space; 117 118 for (i=0;i<am;i++) { 119 apnz = 0; 120 /* diagonal portion of A */ 121 nzi = adi[i+1] - adi[i]; 122 for (j=0; j<nzi; j++){ 123 row = *adj++; 124 pnz = pi_loc[row+1] - pi_loc[row]; 125 Jptr = pj_loc + pi_loc[row]; 126 /* add non-zero cols of P into the sorted linked list lnk */ 127 ierr = PetscLLAdd(pnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr); 128 apnz += nlnk; 129 } 130 /* off-diagonal portion of A */ 131 nzi = aoi[i+1] - aoi[i]; 132 for (j=0; j<nzi; j++){ 133 row = *aoj++; 134 pnz = pi_oth[row+1] - pi_oth[row]; 135 Jptr = pj_oth + pi_oth[row]; 136 ierr = PetscLLAdd(pnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr); 137 apnz += nlnk; 138 } 139 140 api[i+1] = api[i] + apnz; 141 if (ap->abnz_max < apnz) ap->abnz_max = apnz; 142 143 /* if free space is not available, double the total space in the list */ 144 if (current_space->local_remaining<apnz) { 145 ierr = GetMoreSpace(current_space->total_array_size,¤t_space);CHKERRQ(ierr); 146 } 147 148 /* Copy data into free space, then initialize lnk */ 149 ierr = PetscLLClean(pN,pN,apnz,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 150 current_space->array += apnz; 151 current_space->local_used += apnz; 152 current_space->local_remaining -= apnz; 153 } 154 /* Allocate space for apj, initialize apj, and */ 155 /* destroy list of free space and other temporary array(s) */ 156 ierr = PetscMalloc((api[am]+1)*sizeof(PetscInt),&ap->abj);CHKERRQ(ierr); 157 apj = ap->abj; 158 ierr = MakeSpaceContiguous(&free_space,ap->abj);CHKERRQ(ierr); 159 160 /* determine symbolic Co=(p->B)^T*AP - send to others */ 161 /*----------------------------------------------------*/ 162 ierr = MatGetSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr); 163 164 /* then, compute symbolic Co = (p->B)^T*AP */ 165 pon = (p->B)->n; /* total num of rows to be sent to other processors 166 >= (num of nonzero rows of C_seq) - pn */ 167 ierr = PetscMalloc((pon+1)*sizeof(PetscInt),&coi);CHKERRQ(ierr); 168 coi[0] = 0; 169 170 /* set initial free space to be 3*pon*max( nnz(AP) per row) */ 171 nnz = 3*pon*ap->abnz_max + 1; 172 ierr = GetMoreSpace(nnz,&free_space); 173 current_space = free_space; 174 175 for (i=0; i<pon; i++) { 176 nnz = 0; 177 pnz = poti[i+1] - poti[i]; 178 j = pnz; 179 ptJ = potj + poti[i+1]; 180 while (j){/* assume cols are almost in increasing order, starting from its end saves computation */ 181 j--; ptJ--; 182 row = *ptJ; /* row of AP == col of Pot */ 183 apnz = api[row+1] - api[row]; 184 Jptr = apj + api[row]; 185 /* add non-zero cols of AP into the sorted linked list lnk */ 186 ierr = PetscLLAdd(apnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr); 187 nnz += nlnk; 188 } 189 190 /* If free space is not available, double the total space in the list */ 191 if (current_space->local_remaining<nnz) { 192 ierr = GetMoreSpace(current_space->total_array_size,¤t_space);CHKERRQ(ierr); 193 } 194 195 /* Copy data into free space, and zero out denserows */ 196 ierr = PetscLLClean(pN,pN,nnz,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 197 current_space->array += nnz; 198 current_space->local_used += nnz; 199 current_space->local_remaining -= nnz; 200 coi[i+1] = coi[i] + nnz; 201 } 202 ierr = PetscMalloc((coi[pon]+1)*sizeof(PetscInt),&coj);CHKERRQ(ierr); 203 ierr = MakeSpaceContiguous(&free_space,coj);CHKERRQ(ierr); 204 ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr); 205 206 /* send j-array (coj) of Co to other processors */ 207 /*----------------------------------------------*/ 208 /* determine row ownership */ 209 ierr = PetscNew(Mat_Merge_SeqsToMPI,&merge);CHKERRQ(ierr); 210 ierr = PetscMapCreate(comm,&merge->rowmap);CHKERRQ(ierr); 211 ierr = PetscMapSetLocalSize(merge->rowmap,pn);CHKERRQ(ierr); 212 ierr = PetscMapSetType(merge->rowmap,MAP_MPI);CHKERRQ(ierr); 213 ierr = PetscMapGetGlobalRange(merge->rowmap,&owners);CHKERRQ(ierr); 214 215 /* determine the number of messages to send, their lengths */ 216 ierr = PetscMalloc(size*sizeof(PetscMPIInt),&len_si);CHKERRQ(ierr); 217 ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 218 ierr = PetscMalloc(size*sizeof(PetscMPIInt),&merge->len_s);CHKERRQ(ierr); 219 len_s = merge->len_s; 220 merge->nsend = 0; 221 222 ierr = PetscMalloc((size+2)*sizeof(PetscInt),&owners_co);CHKERRQ(ierr); 223 ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 224 225 proc = 0; 226 for (i=0; i<pon; i++){ 227 while (prmap[i] >= owners[proc+1]) proc++; 228 len_si[proc]++; /* num of rows in Co to be sent to [proc] */ 229 len_s[proc] += coi[i+1] - coi[i]; 230 } 231 232 len = 0; /* max length of buf_si[] */ 233 owners_co[0] = 0; 234 for (proc=0; proc<size; proc++){ 235 owners_co[proc+1] = owners_co[proc] + len_si[proc]; 236 if (len_si[proc]){ 237 merge->nsend++; 238 len_si[proc] = 2*(len_si[proc] + 1); 239 len += len_si[proc]; 240 } 241 } 242 243 /* determine the number and length of messages to receive for coi and coj */ 244 ierr = PetscGatherNumberOfMessages(comm,PETSC_NULL,len_s,&merge->nrecv);CHKERRQ(ierr); 245 ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr); 246 247 /* post the Irecv and Isend of coj */ 248 ierr = PetscObjectGetNewTag((PetscObject)merge->rowmap,&tag);CHKERRQ(ierr); 249 ierr = PetscPostIrecvInt(comm,tag,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr); 250 251 ierr = PetscMalloc((merge->nsend+1)*sizeof(MPI_Request),&swaits);CHKERRQ(ierr); 252 253 for (proc=0, k=0; proc<size; proc++){ 254 if (!len_s[proc]) continue; 255 i = owners_co[proc]; 256 ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tag,comm,swaits+k);CHKERRQ(ierr); 257 k++; 258 } 259 260 /* receives and sends of coj are complete */ 261 ierr = PetscMalloc(size*sizeof(MPI_Status),&sstatus);CHKERRQ(ierr); 262 i = merge->nrecv; 263 while (i--) { 264 ierr = MPI_Waitany(merge->nrecv,rwaits,&j,&rstatus);CHKERRQ(ierr); 265 } 266 ierr = PetscFree(rwaits);CHKERRQ(ierr); 267 if (merge->nsend){ 268 ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr); 269 } 270 271 /* send and recv coi */ 272 /*-------------------*/ 273 ierr = PetscPostIrecvInt(comm,tag,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr); 274 275 ierr = PetscMalloc((len+1)*sizeof(PetscInt),&buf_s);CHKERRQ(ierr); 276 buf_si = buf_s; /* points to the beginning of k-th msg to be sent */ 277 for (proc=0,k=0; proc<size; proc++){ 278 if (!len_s[proc]) continue; 279 /* form outgoing message for i-structure: 280 buf_si[0]: nrows to be sent 281 [1:nrows]: row index (global) 282 [nrows+1:2*nrows+1]: i-structure index 283 */ 284 /*-------------------------------------------*/ 285 nrows = len_si[proc]/2 - 1; 286 buf_si_i = buf_si + nrows+1; 287 buf_si[0] = nrows; 288 buf_si_i[0] = 0; 289 nrows = 0; 290 for (i=owners_co[proc]; i<owners_co[proc+1]; i++){ 291 nzi = coi[i+1] - coi[i]; 292 buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */ 293 buf_si[nrows+1] =prmap[i] -owners[proc]; /* local row index */ 294 nrows++; 295 } 296 ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tag,comm,swaits+k);CHKERRQ(ierr); 297 k++; 298 buf_si += len_si[proc]; 299 } 300 i = merge->nrecv; 301 while (i--) { 302 ierr = MPI_Waitany(merge->nrecv,rwaits,&j,&rstatus);CHKERRQ(ierr); 303 } 304 ierr = PetscFree(rwaits);CHKERRQ(ierr); 305 if (merge->nsend){ 306 ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr); 307 } 308 309 ierr = PetscLogInfo(((PetscObject)A,"MatMerge_SeqsToMPI: nsend: %d, nrecv: %d\n",merge->nsend,merge->nrecv));CHKERRQ(ierr); 310 for (i=0; i<merge->nrecv; i++){ 311 ierr = PetscLogInfo(((PetscObject)A,"MatMerge_SeqsToMPI: recv len_ri=%d, len_rj=%d from [%d]\n",len_ri[i],merge->len_r[i],merge->id_r[i]));CHKERRQ(ierr); 312 } 313 314 ierr = PetscFree(len_si);CHKERRQ(ierr); 315 ierr = PetscFree(len_ri);CHKERRQ(ierr); 316 ierr = PetscFree(swaits);CHKERRQ(ierr); 317 ierr = PetscFree(sstatus);CHKERRQ(ierr); 318 ierr = PetscFree(buf_s);CHKERRQ(ierr); 319 320 /* compute the local portion of C (mpi mat) */ 321 /*------------------------------------------*/ 322 ierr = MatGetSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr); 323 324 /* allocate bi array and free space for accumulating nonzero column info */ 325 ierr = PetscMalloc((pn+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr); 326 bi[0] = 0; 327 328 /* set initial free space to be 3*pn*max( nnz(AP) per row) */ 329 nnz = 3*pn*ap->abnz_max + 1; 330 ierr = GetMoreSpace(nnz,&free_space); 331 current_space = free_space; 332 333 ierr = PetscMalloc((3*merge->nrecv+1)*sizeof(PetscInt**),&buf_ri_k);CHKERRQ(ierr); 334 nextrow = buf_ri_k + merge->nrecv; 335 nextci = nextrow + merge->nrecv; 336 for (k=0; k<merge->nrecv; k++){ 337 buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 338 nrows = *buf_ri_k[k]; 339 nextrow[k] = buf_ri_k[k] + 1; /* next row number of k-th recved i-structure */ 340 nextci[k] = buf_ri_k[k] + (nrows + 1);/* poins to the next i-structure of k-th recved i-structure */ 341 } 342 ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr); 343 for (i=0; i<pn; i++) { 344 /* add pdt[i,:]*AP into lnk */ 345 nnz = 0; 346 pnz = pdti[i+1] - pdti[i]; 347 j = pnz; 348 ptJ = pdtj + pdti[i+1]; 349 while (j){/* assume cols are almost in increasing order, starting from its end saves computation */ 350 j--; ptJ--; 351 row = *ptJ; /* row of AP == col of Pt */ 352 apnz = api[row+1] - api[row]; 353 Jptr = apj + api[row]; 354 /* add non-zero cols of AP into the sorted linked list lnk */ 355 ierr = PetscLLAdd(apnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr); 356 nnz += nlnk; 357 } 358 /* add received col data into lnk */ 359 for (k=0; k<merge->nrecv; k++){ /* k-th received message */ 360 if (i == *nextrow[k]) { /* i-th row */ 361 nzi = *(nextci[k]+1) - *nextci[k]; 362 Jptr = buf_rj[k] + *nextci[k]; 363 ierr = PetscLLAdd(nzi,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr); 364 nnz += nlnk; 365 nextrow[k]++; nextci[k]++; 366 } 367 } 368 369 /* if free space is not available, make more free space */ 370 if (current_space->local_remaining<nnz) { 371 ierr = GetMoreSpace(current_space->total_array_size,¤t_space);CHKERRQ(ierr); 372 } 373 /* copy data into free space, then initialize lnk */ 374 ierr = PetscLLClean(pN,pN,nnz,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 375 ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr); 376 current_space->array += nnz; 377 current_space->local_used += nnz; 378 current_space->local_remaining -= nnz; 379 bi[i+1] = bi[i] + nnz; 380 } 381 ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr); 382 ierr = PetscFree(buf_ri_k);CHKERRQ(ierr); 383 384 ierr = PetscMalloc((bi[pn]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr); 385 ierr = MakeSpaceContiguous(&free_space,bj);CHKERRQ(ierr); 386 ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 387 388 /* create symbolic parallel matrix B_mpi */ 389 /*---------------------------------------*/ 390 ierr = MatCreate(comm,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE,&B_mpi);CHKERRQ(ierr); 391 ierr = MatSetType(B_mpi,MATMPIAIJ);CHKERRQ(ierr); 392 ierr = MatMPIAIJSetPreallocation(B_mpi,0,dnz,0,onz);CHKERRQ(ierr); 393 ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 394 395 /* B_mpi is not ready for use - assembly will be done by MatPtAPNumeric() */ 396 B_mpi->assembled = PETSC_FALSE; 397 B_mpi->ops->destroy = MatDestroy_MPIAIJ_SeqsToMPI; 398 merge->bi = bi; 399 merge->bj = bj; 400 merge->coi = coi; 401 merge->coj = coj; 402 merge->buf_ri = buf_ri; 403 merge->buf_rj = buf_rj; 404 merge->owners_co = owners_co; 405 406 /* attach the supporting struct to B_mpi for reuse */ 407 ierr = PetscObjectContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr); 408 ierr = PetscObjectContainerSetPointer(container,merge);CHKERRQ(ierr); 409 ierr = PetscObjectCompose((PetscObject)B_mpi,"MatMergeSeqsToMPI",(PetscObject)container);CHKERRQ(ierr); 410 *C = B_mpi; 411 PetscFunctionReturn(0); 412 } 413 414 #undef __FUNCT__ 415 #define __FUNCT__ "MatPtAPNumeric_MPIAIJ_MPIAIJ" 416 PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C) 417 { 418 PetscErrorCode ierr; 419 Mat_Merge_SeqsToMPI *merge; 420 Mat_MatMatMultMPI *ap; 421 PetscObjectContainer cont_merge,cont_ptap; 422 Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data; 423 Mat_SeqAIJ *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data; 424 Mat_SeqAIJ *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data; 425 Mat_SeqAIJ *p_loc,*p_oth; 426 PetscInt *adi=ad->i,*aoi=ao->i,*adj=ad->j,*aoj=ao->j,*apJ,nextp,flops=0; 427 PetscInt *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pJ,*pj; 428 PetscInt i,j,k,anz,pnz,apnz,nextap,row,*cj; 429 MatScalar *ada=ad->a,*aoa=ao->a,*apa,*pa,*ca,*pa_loc,*pa_oth; 430 PetscInt am=A->m,cm=C->m,pon=(p->B)->n; 431 MPI_Comm comm=C->comm; 432 PetscMPIInt size,rank,taga,*len_s; 433 PetscInt *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci; 434 PetscInt **buf_ri,**buf_rj; 435 PetscInt cnz=0,*bj_i,*bi,*bj,bnz,nextcj; /* bi,bj,ba: local array of C(mpi mat) */ 436 MPI_Request *s_waits,*r_waits; 437 MPI_Status *status; 438 MatScalar **abuf_r,*ba_i,*pA,*coa,*ba; 439 PetscInt *api,*apj,*coi,*coj; 440 PetscInt *poJ=po->j,*pdJ=pd->j,pcstart=p->cstart,pcend=p->cend; 441 442 PetscFunctionBegin; 443 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 444 ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 445 446 ierr = PetscObjectQuery((PetscObject)C,"MatMergeSeqsToMPI",(PetscObject *)&cont_merge);CHKERRQ(ierr); 447 if (cont_merge) { 448 ierr = PetscObjectContainerGetPointer(cont_merge,(void **)&merge);CHKERRQ(ierr); 449 } else { 450 SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Matrix C does not posses an object container"); 451 } 452 453 ierr = PetscObjectQuery((PetscObject)P,"Mat_MatMatMultMPI",(PetscObject *)&cont_ptap);CHKERRQ(ierr); 454 if (cont_ptap) { 455 ierr = PetscObjectContainerGetPointer(cont_ptap,(void **)&ap);CHKERRQ(ierr); 456 if (ap->reuse == MAT_INITIAL_MATRIX){ 457 ap->reuse = MAT_REUSE_MATRIX; 458 } else { /* update numerical values of P_oth and P_loc */ 459 ierr = MatGetBrowsOfAoCols(A,P,MAT_REUSE_MATRIX,&ap->startsj,&ap->bufa,&ap->B_oth);CHKERRQ(ierr); 460 ierr = MatGetLocalMat(P,MAT_REUSE_MATRIX,&ap->B_loc);CHKERRQ(ierr); 461 } 462 } else { 463 SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Matrix P does not posses an object container"); 464 } 465 466 /* get data from symbolic products */ 467 p_loc = (Mat_SeqAIJ*)(ap->B_loc)->data; 468 p_oth = (Mat_SeqAIJ*)(ap->B_oth)->data; 469 pi_loc=p_loc->i; pj_loc=p_loc->j; pJ=pj_loc; pa_loc=p_loc->a,pA=pa_loc; 470 pi_oth=p_oth->i; pj_oth=p_oth->j; pa_oth=p_oth->a; 471 472 coi = merge->coi; coj = merge->coj; 473 ierr = PetscMalloc((coi[pon]+1)*sizeof(MatScalar),&coa);CHKERRQ(ierr); 474 ierr = PetscMemzero(coa,coi[pon]*sizeof(MatScalar));CHKERRQ(ierr); 475 476 bi = merge->bi; bj = merge->bj; 477 ierr = PetscMapGetGlobalRange(merge->rowmap,&owners);CHKERRQ(ierr); 478 ierr = PetscMalloc((bi[cm]+1)*sizeof(MatScalar),&ba);CHKERRQ(ierr); 479 ierr = PetscMemzero(ba,bi[cm]*sizeof(MatScalar));CHKERRQ(ierr); 480 481 /* get data from symbolic A*P */ 482 ierr = PetscMalloc((ap->abnz_max+1)*sizeof(MatScalar),&apa);CHKERRQ(ierr); 483 ierr = PetscMemzero(apa,ap->abnz_max*sizeof(MatScalar));CHKERRQ(ierr); 484 485 /* compute numeric C_seq=P_loc^T*A_loc*P */ 486 api = ap->abi; apj = ap->abj; 487 for (i=0;i<am;i++) { 488 /* form i-th sparse row of A*P */ 489 apnz = api[i+1] - api[i]; 490 apJ = apj + api[i]; 491 /* diagonal portion of A */ 492 anz = adi[i+1] - adi[i]; 493 for (j=0;j<anz;j++) { 494 row = *adj++; 495 pnz = pi_loc[row+1] - pi_loc[row]; 496 pj = pj_loc + pi_loc[row]; 497 pa = pa_loc + pi_loc[row]; 498 nextp = 0; 499 for (k=0; nextp<pnz; k++) { 500 if (apJ[k] == pj[nextp]) { /* col of AP == col of P */ 501 apa[k] += (*ada)*pa[nextp++]; 502 } 503 } 504 flops += 2*pnz; 505 ada++; 506 } 507 /* off-diagonal portion of A */ 508 anz = aoi[i+1] - aoi[i]; 509 for (j=0; j<anz; j++) { 510 row = *aoj++; 511 pnz = pi_oth[row+1] - pi_oth[row]; 512 pj = pj_oth + pi_oth[row]; 513 pa = pa_oth + pi_oth[row]; 514 nextp = 0; 515 for (k=0; nextp<pnz; k++) { 516 if (apJ[k] == pj[nextp]) { /* col of AP == col of P */ 517 apa[k] += (*aoa)*pa[nextp++]; 518 } 519 } 520 flops += 2*pnz; 521 aoa++; 522 } 523 524 /* Compute P_loc[i,:]^T*AP[i,:] using outer product */ 525 pnz = pi_loc[i+1] - pi_loc[i]; 526 for (j=0; j<pnz; j++) { 527 nextap = 0; 528 row = *pJ++; /* global index */ 529 if (row < pcstart || row >=pcend) { /* put the value into Co */ 530 cj = coj + coi[*poJ]; 531 ca = coa + coi[*poJ++]; 532 } else { /* put the value into Cd */ 533 cj = bj + bi[*pdJ]; 534 ca = ba + bi[*pdJ++]; 535 } 536 for (k=0; nextap<apnz; k++) { 537 if (cj[k]==apJ[nextap]) ca[k] += (*pA)*apa[nextap++]; 538 } 539 flops += 2*apnz; 540 pA++; 541 } 542 543 /* zero the current row info for A*P */ 544 ierr = PetscMemzero(apa,apnz*sizeof(MatScalar));CHKERRQ(ierr); 545 } 546 ierr = PetscFree(apa);CHKERRQ(ierr); 547 548 /* send and recv matrix values */ 549 /*-----------------------------*/ 550 buf_ri = merge->buf_ri; 551 buf_rj = merge->buf_rj; 552 len_s = merge->len_s; 553 ierr = PetscObjectGetNewTag((PetscObject)merge->rowmap,&taga);CHKERRQ(ierr); 554 ierr = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr); 555 556 ierr = PetscMalloc((merge->nsend+1)*sizeof(MPI_Request),&s_waits);CHKERRQ(ierr); 557 for (proc=0,k=0; proc<size; proc++){ 558 if (!len_s[proc]) continue; 559 i = merge->owners_co[proc]; 560 ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr); 561 k++; 562 } 563 ierr = PetscMalloc(size*sizeof(MPI_Status),&status);CHKERRQ(ierr); 564 ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr); 565 ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr); 566 ierr = PetscFree(status);CHKERRQ(ierr); 567 568 ierr = PetscFree(s_waits);CHKERRQ(ierr); 569 ierr = PetscFree(r_waits);CHKERRQ(ierr); 570 ierr = PetscFree(coa);CHKERRQ(ierr); 571 572 /* insert local and received values into C */ 573 /*-----------------------------------------*/ 574 ierr = PetscMalloc((3*merge->nrecv+1)*sizeof(PetscInt**),&buf_ri_k);CHKERRQ(ierr); 575 nextrow = buf_ri_k + merge->nrecv; 576 nextci = nextrow + merge->nrecv; 577 578 for (k=0; k<merge->nrecv; k++){ 579 buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 580 nrows = *(buf_ri_k[k]); 581 nextrow[k] = buf_ri_k[k]+1; /* next row number of k-th recved i-structure */ 582 nextci[k] = buf_ri_k[k] + (nrows + 1);/* poins to the next i-structure of k-th recved i-structure */ 583 } 584 585 for (i=0; i<cm; i++) { 586 row = owners[rank] + i; /* global row index of C_seq */ 587 bj_i = bj + bi[i]; /* col indices of the i-th row of C */ 588 ba_i = ba + bi[i]; 589 bnz = bi[i+1] - bi[i]; 590 /* add received vals into ba */ 591 for (k=0; k<merge->nrecv; k++){ /* k-th received message */ 592 /* i-th row */ 593 if (i == *nextrow[k]) { 594 cnz = *(nextci[k]+1) - *nextci[k]; 595 cj = buf_rj[k] + *(nextci[k]); 596 ca = abuf_r[k] + *(nextci[k]); 597 nextcj = 0; 598 for (j=0; nextcj<cnz; j++){ 599 if (bj_i[j] == cj[nextcj]){ /* bcol == ccol */ 600 ba_i[j] += ca[nextcj++]; 601 } 602 } 603 nextrow[k]++; nextci[k]++; 604 } 605 } 606 ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr); 607 flops += 2*cnz; 608 } 609 ierr = MatSetOption(C,MAT_COLUMNS_SORTED);CHKERRQ(ierr); /* -cause delay? */ 610 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 611 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 612 613 ierr = PetscFree(ba);CHKERRQ(ierr); 614 ierr = PetscFree(abuf_r);CHKERRQ(ierr); 615 ierr = PetscFree(buf_ri_k);CHKERRQ(ierr); 616 ierr = PetscLogFlops(flops);CHKERRQ(ierr); 617 PetscFunctionReturn(0); 618 } 619