1 /* 2 Routines to compute overlapping regions of a parallel MPI matrix 3 and to find submatrices that were shared across processors. 4 */ 5 #include "src/mat/impls/aij/mpi/mpiaij.h" 6 #include "petscbt.h" 7 8 static PetscErrorCode MatIncreaseOverlap_MPIAIJ_Once(Mat,int,IS *); 9 static PetscErrorCode MatIncreaseOverlap_MPIAIJ_Local(Mat,int,char **,int*,int**); 10 static PetscErrorCode MatIncreaseOverlap_MPIAIJ_Receive(Mat,int,int **,int**,int*); 11 EXTERN PetscErrorCode MatGetRow_MPIAIJ(Mat,int,int*,int**,PetscScalar**); 12 EXTERN PetscErrorCode MatRestoreRow_MPIAIJ(Mat,int,int*,int**,PetscScalar**); 13 14 #undef __FUNCT__ 15 #define __FUNCT__ "MatIncreaseOverlap_MPIAIJ" 16 PetscErrorCode MatIncreaseOverlap_MPIAIJ(Mat C,int imax,IS is[],int ov) 17 { 18 PetscErrorCode ierr; 19 int i; 20 21 PetscFunctionBegin; 22 if (ov < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap specified"); 23 for (i=0; i<ov; ++i) { 24 ierr = MatIncreaseOverlap_MPIAIJ_Once(C,imax,is);CHKERRQ(ierr); 25 } 26 PetscFunctionReturn(0); 27 } 28 29 /* 30 Sample message format: 31 If a processor A wants processor B to process some elements corresponding 32 to index sets is[1],is[5] 33 mesg [0] = 2 (no of index sets in the mesg) 34 ----------- 35 mesg [1] = 1 => is[1] 36 mesg [2] = sizeof(is[1]); 37 ----------- 38 mesg [3] = 5 => is[5] 39 mesg [4] = sizeof(is[5]); 40 ----------- 41 mesg [5] 42 mesg [n] datas[1] 43 ----------- 44 mesg[n+1] 45 mesg[m] data(is[5]) 46 ----------- 47 48 Notes: 49 nrqs - no of requests sent (or to be sent out) 50 nrqr - no of requests recieved (which have to be or which have been processed 51 */ 52 #undef __FUNCT__ 53 #define __FUNCT__ "MatIncreaseOverlap_MPIAIJ_Once" 54 static PetscErrorCode MatIncreaseOverlap_MPIAIJ_Once(Mat C,int imax,IS is[]) 55 { 56 Mat_MPIAIJ *c = (Mat_MPIAIJ*)C->data; 57 int **idx,*n,*w1,*w2,*w3,*w4,*rtable,**data,len,*idx_i; 58 PetscErrorCode ierr; 59 int size,rank,m,i,j,k,**rbuf,row,proc,nrqs,msz,**outdat,**ptr; 60 int *ctr,*pa,*tmp,nrqr,*isz,*isz1,**xdata,**rbuf2; 61 int *onodes1,*olengths1,tag1,tag2,*onodes2,*olengths2; 62 PetscBT *table; 63 MPI_Comm comm; 64 MPI_Request *s_waits1,*r_waits1,*s_waits2,*r_waits2; 65 MPI_Status *s_status,*recv_status; 66 67 PetscFunctionBegin; 68 comm = C->comm; 69 size = c->size; 70 rank = c->rank; 71 m = C->M; 72 73 ierr = PetscObjectGetNewTag((PetscObject)C,&tag1);CHKERRQ(ierr); 74 ierr = PetscObjectGetNewTag((PetscObject)C,&tag2);CHKERRQ(ierr); 75 76 len = (imax+1)*sizeof(int*)+ (imax + m)*sizeof(int); 77 ierr = PetscMalloc(len,&idx);CHKERRQ(ierr); 78 n = (int*)(idx + imax); 79 rtable = n + imax; 80 81 for (i=0; i<imax; i++) { 82 ierr = ISGetIndices(is[i],&idx[i]);CHKERRQ(ierr); 83 ierr = ISGetLocalSize(is[i],&n[i]);CHKERRQ(ierr); 84 } 85 86 /* Create hash table for the mapping :row -> proc*/ 87 for (i=0,j=0; i<size; i++) { 88 len = c->rowners[i+1]; 89 for (; j<len; j++) { 90 rtable[j] = i; 91 } 92 } 93 94 /* evaluate communication - mesg to who,length of mesg, and buffer space 95 required. Based on this, buffers are allocated, and data copied into them*/ 96 ierr = PetscMalloc(size*4*sizeof(int),&w1);CHKERRQ(ierr);/* mesg size */ 97 w2 = w1 + size; /* if w2[i] marked, then a message to proc i*/ 98 w3 = w2 + size; /* no of IS that needs to be sent to proc i */ 99 w4 = w3 + size; /* temp work space used in determining w1, w2, w3 */ 100 ierr = PetscMemzero(w1,size*3*sizeof(int));CHKERRQ(ierr); /* initialise work vector*/ 101 for (i=0; i<imax; i++) { 102 ierr = PetscMemzero(w4,size*sizeof(int));CHKERRQ(ierr); /* initialise work vector*/ 103 idx_i = idx[i]; 104 len = n[i]; 105 for (j=0; j<len; j++) { 106 row = idx_i[j]; 107 if (row < 0) { 108 SETERRQ(1,"Index set cannot have negative entries"); 109 } 110 proc = rtable[row]; 111 w4[proc]++; 112 } 113 for (j=0; j<size; j++){ 114 if (w4[j]) { w1[j] += w4[j]; w3[j]++;} 115 } 116 } 117 118 nrqs = 0; /* no of outgoing messages */ 119 msz = 0; /* total mesg length (for all proc */ 120 w1[rank] = 0; /* no mesg sent to intself */ 121 w3[rank] = 0; 122 for (i=0; i<size; i++) { 123 if (w1[i]) {w2[i] = 1; nrqs++;} /* there exists a message to proc i */ 124 } 125 /* pa - is list of processors to communicate with */ 126 ierr = PetscMalloc((nrqs+1)*sizeof(int),&pa);CHKERRQ(ierr); 127 for (i=0,j=0; i<size; i++) { 128 if (w1[i]) {pa[j] = i; j++;} 129 } 130 131 /* Each message would have a header = 1 + 2*(no of IS) + data */ 132 for (i=0; i<nrqs; i++) { 133 j = pa[i]; 134 w1[j] += w2[j] + 2*w3[j]; 135 msz += w1[j]; 136 } 137 138 /* Determine the number of messages to expect, their lengths, from from-ids */ 139 ierr = PetscGatherNumberOfMessages(comm,w2,w1,&nrqr);CHKERRQ(ierr); 140 ierr = PetscGatherMessageLengths(comm,nrqs,nrqr,w1,&onodes1,&olengths1);CHKERRQ(ierr); 141 142 /* Now post the Irecvs corresponding to these messages */ 143 ierr = PetscPostIrecvInt(comm,tag1,nrqr,onodes1,olengths1,&rbuf,&r_waits1);CHKERRQ(ierr); 144 145 /* Allocate Memory for outgoing messages */ 146 len = 2*size*sizeof(int*) + (size+msz)*sizeof(int); 147 ierr = PetscMalloc(len,&outdat);CHKERRQ(ierr); 148 ptr = outdat + size; /* Pointers to the data in outgoing buffers */ 149 ierr = PetscMemzero(outdat,2*size*sizeof(int*));CHKERRQ(ierr); 150 tmp = (int*)(outdat + 2*size); 151 ctr = tmp + msz; 152 153 { 154 int *iptr = tmp,ict = 0; 155 for (i=0; i<nrqs; i++) { 156 j = pa[i]; 157 iptr += ict; 158 outdat[j] = iptr; 159 ict = w1[j]; 160 } 161 } 162 163 /* Form the outgoing messages */ 164 /*plug in the headers*/ 165 for (i=0; i<nrqs; i++) { 166 j = pa[i]; 167 outdat[j][0] = 0; 168 ierr = PetscMemzero(outdat[j]+1,2*w3[j]*sizeof(int));CHKERRQ(ierr); 169 ptr[j] = outdat[j] + 2*w3[j] + 1; 170 } 171 172 /* Memory for doing local proc's work*/ 173 { 174 int *d_p; 175 char *t_p; 176 177 len = (imax)*(sizeof(PetscBT) + sizeof(int*)+ sizeof(int)) + 178 (m)*imax*sizeof(int) + (m/PETSC_BITS_PER_BYTE+1)*imax*sizeof(char) + 1; 179 ierr = PetscMalloc(len,&table);CHKERRQ(ierr); 180 ierr = PetscMemzero(table,len);CHKERRQ(ierr); 181 data = (int **)(table + imax); 182 isz = (int *)(data + imax); 183 d_p = (int *)(isz + imax); 184 t_p = (char *)(d_p + m*imax); 185 for (i=0; i<imax; i++) { 186 table[i] = t_p + (m/PETSC_BITS_PER_BYTE+1)*i; 187 data[i] = d_p + (m)*i; 188 } 189 } 190 191 /* Parse the IS and update local tables and the outgoing buf with the data*/ 192 { 193 int n_i,*data_i,isz_i,*outdat_j,ctr_j; 194 PetscBT table_i; 195 196 for (i=0; i<imax; i++) { 197 ierr = PetscMemzero(ctr,size*sizeof(int));CHKERRQ(ierr); 198 n_i = n[i]; 199 table_i = table[i]; 200 idx_i = idx[i]; 201 data_i = data[i]; 202 isz_i = isz[i]; 203 for (j=0; j<n_i; j++) { /* parse the indices of each IS */ 204 row = idx_i[j]; 205 proc = rtable[row]; 206 if (proc != rank) { /* copy to the outgoing buffer */ 207 ctr[proc]++; 208 *ptr[proc] = row; 209 ptr[proc]++; 210 } else { /* Update the local table */ 211 if (!PetscBTLookupSet(table_i,row)) { data_i[isz_i++] = row;} 212 } 213 } 214 /* Update the headers for the current IS */ 215 for (j=0; j<size; j++) { /* Can Optimise this loop by using pa[] */ 216 if ((ctr_j = ctr[j])) { 217 outdat_j = outdat[j]; 218 k = ++outdat_j[0]; 219 outdat_j[2*k] = ctr_j; 220 outdat_j[2*k-1] = i; 221 } 222 } 223 isz[i] = isz_i; 224 } 225 } 226 227 228 229 /* Now post the sends */ 230 ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Request),&s_waits1);CHKERRQ(ierr); 231 for (i=0; i<nrqs; ++i) { 232 j = pa[i]; 233 ierr = MPI_Isend(outdat[j],w1[j],MPI_INT,j,tag1,comm,s_waits1+i);CHKERRQ(ierr); 234 } 235 236 /* No longer need the original indices*/ 237 for (i=0; i<imax; ++i) { 238 ierr = ISRestoreIndices(is[i],idx+i);CHKERRQ(ierr); 239 } 240 ierr = PetscFree(idx);CHKERRQ(ierr); 241 242 for (i=0; i<imax; ++i) { 243 ierr = ISDestroy(is[i]);CHKERRQ(ierr); 244 } 245 246 /* Do Local work*/ 247 ierr = MatIncreaseOverlap_MPIAIJ_Local(C,imax,table,isz,data);CHKERRQ(ierr); 248 249 /* Receive messages*/ 250 ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Status),&recv_status);CHKERRQ(ierr); 251 ierr = MPI_Waitall(nrqr,r_waits1,recv_status);CHKERRQ(ierr); 252 253 ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Status),&s_status);CHKERRQ(ierr); 254 ierr = MPI_Waitall(nrqs,s_waits1,s_status);CHKERRQ(ierr); 255 256 /* Phase 1 sends are complete - deallocate buffers */ 257 ierr = PetscFree(outdat);CHKERRQ(ierr); 258 ierr = PetscFree(w1);CHKERRQ(ierr); 259 260 ierr = PetscMalloc((nrqr+1)*sizeof(int*),&xdata);CHKERRQ(ierr); 261 ierr = PetscMalloc((nrqr+1)*sizeof(int),&isz1);CHKERRQ(ierr); 262 ierr = MatIncreaseOverlap_MPIAIJ_Receive(C,nrqr,rbuf,xdata,isz1);CHKERRQ(ierr); 263 ierr = PetscFree(rbuf);CHKERRQ(ierr); 264 265 266 /* Send the data back*/ 267 /* Do a global reduction to know the buffer space req for incoming messages*/ 268 { 269 int *rw1; 270 271 ierr = PetscMalloc(size*sizeof(int),&rw1);CHKERRQ(ierr); 272 ierr = PetscMemzero(rw1,size*sizeof(int));CHKERRQ(ierr); 273 274 for (i=0; i<nrqr; ++i) { 275 proc = recv_status[i].MPI_SOURCE; 276 if (proc != onodes1[i]) SETERRQ(1,"MPI_SOURCE mismatch"); 277 rw1[proc] = isz1[i]; 278 } 279 ierr = PetscFree(onodes1);CHKERRQ(ierr); 280 ierr = PetscFree(olengths1);CHKERRQ(ierr); 281 282 /* Determine the number of messages to expect, their lengths, from from-ids */ 283 ierr = PetscGatherMessageLengths(comm,nrqr,nrqs,rw1,&onodes2,&olengths2);CHKERRQ(ierr); 284 PetscFree(rw1); 285 } 286 /* Now post the Irecvs corresponding to these messages */ 287 ierr = PetscPostIrecvInt(comm,tag2,nrqs,onodes2,olengths2,&rbuf2,&r_waits2);CHKERRQ(ierr); 288 289 /* Now post the sends */ 290 ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Request),&s_waits2);CHKERRQ(ierr); 291 for (i=0; i<nrqr; ++i) { 292 j = recv_status[i].MPI_SOURCE; 293 ierr = MPI_Isend(xdata[i],isz1[i],MPI_INT,j,tag2,comm,s_waits2+i);CHKERRQ(ierr); 294 } 295 296 /* receive work done on other processors*/ 297 { 298 int idex,is_no,ct1,max,*rbuf2_i,isz_i,*data_i,jmax; 299 PetscBT table_i; 300 MPI_Status *status2; 301 302 ierr = PetscMalloc((PetscMax(nrqr,nrqs)+1)*sizeof(MPI_Status),&status2);CHKERRQ(ierr); 303 for (i=0; i<nrqs; ++i) { 304 ierr = MPI_Waitany(nrqs,r_waits2,&idex,status2+i);CHKERRQ(ierr); 305 /* Process the message*/ 306 rbuf2_i = rbuf2[idex]; 307 ct1 = 2*rbuf2_i[0]+1; 308 jmax = rbuf2[idex][0]; 309 for (j=1; j<=jmax; j++) { 310 max = rbuf2_i[2*j]; 311 is_no = rbuf2_i[2*j-1]; 312 isz_i = isz[is_no]; 313 data_i = data[is_no]; 314 table_i = table[is_no]; 315 for (k=0; k<max; k++,ct1++) { 316 row = rbuf2_i[ct1]; 317 if (!PetscBTLookupSet(table_i,row)) { data_i[isz_i++] = row;} 318 } 319 isz[is_no] = isz_i; 320 } 321 } 322 323 ierr = MPI_Waitall(nrqr,s_waits2,status2);CHKERRQ(ierr); 324 ierr = PetscFree(status2);CHKERRQ(ierr); 325 } 326 327 for (i=0; i<imax; ++i) { 328 ierr = ISCreateGeneral(PETSC_COMM_SELF,isz[i],data[i],is+i);CHKERRQ(ierr); 329 } 330 331 ierr = PetscFree(onodes2);CHKERRQ(ierr); 332 ierr = PetscFree(olengths2);CHKERRQ(ierr); 333 334 ierr = PetscFree(pa);CHKERRQ(ierr); 335 ierr = PetscFree(rbuf2);CHKERRQ(ierr); 336 ierr = PetscFree(s_waits1);CHKERRQ(ierr); 337 ierr = PetscFree(r_waits1);CHKERRQ(ierr); 338 ierr = PetscFree(s_waits2);CHKERRQ(ierr); 339 ierr = PetscFree(r_waits2);CHKERRQ(ierr); 340 ierr = PetscFree(table);CHKERRQ(ierr); 341 ierr = PetscFree(s_status);CHKERRQ(ierr); 342 ierr = PetscFree(recv_status);CHKERRQ(ierr); 343 ierr = PetscFree(xdata[0]);CHKERRQ(ierr); 344 ierr = PetscFree(xdata);CHKERRQ(ierr); 345 ierr = PetscFree(isz1);CHKERRQ(ierr); 346 PetscFunctionReturn(0); 347 } 348 349 #undef __FUNCT__ 350 #define __FUNCT__ "MatIncreaseOverlap_MPIAIJ_Local" 351 /* 352 MatIncreaseOverlap_MPIAIJ_Local - Called by MatincreaseOverlap, to do 353 the work on the local processor. 354 355 Inputs: 356 C - MAT_MPIAIJ; 357 imax - total no of index sets processed at a time; 358 table - an array of char - size = m bits. 359 360 Output: 361 isz - array containing the count of the solution elements correspondign 362 to each index set; 363 data - pointer to the solutions 364 */ 365 static PetscErrorCode MatIncreaseOverlap_MPIAIJ_Local(Mat C,int imax,PetscBT *table,int *isz,int **data) 366 { 367 Mat_MPIAIJ *c = (Mat_MPIAIJ*)C->data; 368 Mat A = c->A,B = c->B; 369 Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*b = (Mat_SeqAIJ*)B->data; 370 int start,end,val,max,rstart,cstart,*ai,*aj; 371 int *bi,*bj,*garray,i,j,k,row,*data_i,isz_i; 372 PetscBT table_i; 373 374 PetscFunctionBegin; 375 rstart = c->rstart; 376 cstart = c->cstart; 377 ai = a->i; 378 aj = a->j; 379 bi = b->i; 380 bj = b->j; 381 garray = c->garray; 382 383 384 for (i=0; i<imax; i++) { 385 data_i = data[i]; 386 table_i = table[i]; 387 isz_i = isz[i]; 388 for (j=0,max=isz[i]; j<max; j++) { 389 row = data_i[j] - rstart; 390 start = ai[row]; 391 end = ai[row+1]; 392 for (k=start; k<end; k++) { /* Amat */ 393 val = aj[k] + cstart; 394 if (!PetscBTLookupSet(table_i,val)) { data_i[isz_i++] = val;} 395 } 396 start = bi[row]; 397 end = bi[row+1]; 398 for (k=start; k<end; k++) { /* Bmat */ 399 val = garray[bj[k]]; 400 if (!PetscBTLookupSet(table_i,val)) { data_i[isz_i++] = val;} 401 } 402 } 403 isz[i] = isz_i; 404 } 405 PetscFunctionReturn(0); 406 } 407 408 #undef __FUNCT__ 409 #define __FUNCT__ "MatIncreaseOverlap_MPIAIJ_Receive" 410 /* 411 MatIncreaseOverlap_MPIAIJ_Receive - Process the recieved messages, 412 and return the output 413 414 Input: 415 C - the matrix 416 nrqr - no of messages being processed. 417 rbuf - an array of pointers to the recieved requests 418 419 Output: 420 xdata - array of messages to be sent back 421 isz1 - size of each message 422 423 For better efficiency perhaps we should malloc seperately each xdata[i], 424 then if a remalloc is required we need only copy the data for that one row 425 rather then all previous rows as it is now where a single large chunck of 426 memory is used. 427 428 */ 429 static PetscErrorCode MatIncreaseOverlap_MPIAIJ_Receive(Mat C,int nrqr,int **rbuf,int **xdata,int * isz1) 430 { 431 Mat_MPIAIJ *c = (Mat_MPIAIJ*)C->data; 432 Mat A = c->A,B = c->B; 433 Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*b = (Mat_SeqAIJ*)B->data; 434 PetscErrorCode ierr; 435 int rstart,cstart,*ai,*aj,*bi,*bj,*garray,i,j,k; 436 int row,total_sz,ct,ct1,ct2,ct3,mem_estimate,oct2,l,start,end; 437 int val,max1,max2,rank,m,no_malloc =0,*tmp,new_estimate,ctr; 438 int *rbuf_i,kmax,rbuf_0; 439 PetscBT xtable; 440 441 PetscFunctionBegin; 442 rank = c->rank; 443 m = C->M; 444 rstart = c->rstart; 445 cstart = c->cstart; 446 ai = a->i; 447 aj = a->j; 448 bi = b->i; 449 bj = b->j; 450 garray = c->garray; 451 452 453 for (i=0,ct=0,total_sz=0; i<nrqr; ++i) { 454 rbuf_i = rbuf[i]; 455 rbuf_0 = rbuf_i[0]; 456 ct += rbuf_0; 457 for (j=1; j<=rbuf_0; j++) { total_sz += rbuf_i[2*j]; } 458 } 459 460 if (C->m) max1 = ct*(a->nz + b->nz)/C->m; 461 else max1 = 1; 462 mem_estimate = 3*((total_sz > max1 ? total_sz : max1)+1); 463 ierr = PetscMalloc(mem_estimate*sizeof(int),&xdata[0]);CHKERRQ(ierr); 464 ++no_malloc; 465 ierr = PetscBTCreate(m,xtable);CHKERRQ(ierr); 466 ierr = PetscMemzero(isz1,nrqr*sizeof(int));CHKERRQ(ierr); 467 468 ct3 = 0; 469 for (i=0; i<nrqr; i++) { /* for easch mesg from proc i */ 470 rbuf_i = rbuf[i]; 471 rbuf_0 = rbuf_i[0]; 472 ct1 = 2*rbuf_0+1; 473 ct2 = ct1; 474 ct3 += ct1; 475 for (j=1; j<=rbuf_0; j++) { /* for each IS from proc i*/ 476 ierr = PetscBTMemzero(m,xtable);CHKERRQ(ierr); 477 oct2 = ct2; 478 kmax = rbuf_i[2*j]; 479 for (k=0; k<kmax; k++,ct1++) { 480 row = rbuf_i[ct1]; 481 if (!PetscBTLookupSet(xtable,row)) { 482 if (!(ct3 < mem_estimate)) { 483 new_estimate = (int)(1.5*mem_estimate)+1; 484 ierr = PetscMalloc(new_estimate*sizeof(int),&tmp);CHKERRQ(ierr); 485 ierr = PetscMemcpy(tmp,xdata[0],mem_estimate*sizeof(int));CHKERRQ(ierr); 486 ierr = PetscFree(xdata[0]);CHKERRQ(ierr); 487 xdata[0] = tmp; 488 mem_estimate = new_estimate; ++no_malloc; 489 for (ctr=1; ctr<=i; ctr++) { xdata[ctr] = xdata[ctr-1] + isz1[ctr-1];} 490 } 491 xdata[i][ct2++] = row; 492 ct3++; 493 } 494 } 495 for (k=oct2,max2=ct2; k<max2; k++) { 496 row = xdata[i][k] - rstart; 497 start = ai[row]; 498 end = ai[row+1]; 499 for (l=start; l<end; l++) { 500 val = aj[l] + cstart; 501 if (!PetscBTLookupSet(xtable,val)) { 502 if (!(ct3 < mem_estimate)) { 503 new_estimate = (int)(1.5*mem_estimate)+1; 504 ierr = PetscMalloc(new_estimate*sizeof(int),&tmp);CHKERRQ(ierr); 505 ierr = PetscMemcpy(tmp,xdata[0],mem_estimate*sizeof(int));CHKERRQ(ierr); 506 ierr = PetscFree(xdata[0]);CHKERRQ(ierr); 507 xdata[0] = tmp; 508 mem_estimate = new_estimate; ++no_malloc; 509 for (ctr=1; ctr<=i; ctr++) { xdata[ctr] = xdata[ctr-1] + isz1[ctr-1];} 510 } 511 xdata[i][ct2++] = val; 512 ct3++; 513 } 514 } 515 start = bi[row]; 516 end = bi[row+1]; 517 for (l=start; l<end; l++) { 518 val = garray[bj[l]]; 519 if (!PetscBTLookupSet(xtable,val)) { 520 if (!(ct3 < mem_estimate)) { 521 new_estimate = (int)(1.5*mem_estimate)+1; 522 ierr = PetscMalloc(new_estimate*sizeof(int),&tmp);CHKERRQ(ierr); 523 ierr = PetscMemcpy(tmp,xdata[0],mem_estimate*sizeof(int));CHKERRQ(ierr); 524 ierr = PetscFree(xdata[0]);CHKERRQ(ierr); 525 xdata[0] = tmp; 526 mem_estimate = new_estimate; ++no_malloc; 527 for (ctr =1; ctr <=i; ctr++) { xdata[ctr] = xdata[ctr-1] + isz1[ctr-1];} 528 } 529 xdata[i][ct2++] = val; 530 ct3++; 531 } 532 } 533 } 534 /* Update the header*/ 535 xdata[i][2*j] = ct2 - oct2; /* Undo the vector isz1 and use only a var*/ 536 xdata[i][2*j-1] = rbuf_i[2*j-1]; 537 } 538 xdata[i][0] = rbuf_0; 539 xdata[i+1] = xdata[i] + ct2; 540 isz1[i] = ct2; /* size of each message */ 541 } 542 ierr = PetscBTDestroy(xtable);CHKERRQ(ierr); 543 PetscLogInfo(0,"MatIncreaseOverlap_MPIAIJ:[%d] Allocated %d bytes, required %d bytes, no of mallocs = %d\n",rank,mem_estimate, ct3,no_malloc); 544 PetscFunctionReturn(0); 545 } 546 /* -------------------------------------------------------------------------*/ 547 EXTERN PetscErrorCode MatGetSubMatrices_MPIAIJ_Local(Mat,int,const IS[],const IS[],MatReuse,Mat*); 548 EXTERN PetscErrorCode MatAssemblyEnd_SeqAIJ(Mat,MatAssemblyType); 549 /* 550 Every processor gets the entire matrix 551 */ 552 #undef __FUNCT__ 553 #define __FUNCT__ "MatGetSubMatrix_MPIAIJ_All" 554 PetscErrorCode MatGetSubMatrix_MPIAIJ_All(Mat A,MatReuse scall,Mat *Bin[]) 555 { 556 Mat B; 557 Mat_MPIAIJ *a = (Mat_MPIAIJ *)A->data; 558 Mat_SeqAIJ *b,*ad = (Mat_SeqAIJ*)a->A->data,*bd = (Mat_SeqAIJ*)a->B->data; 559 PetscErrorCode ierr; 560 int sendcount,*recvcounts = 0,*displs = 0,size,i,*rstarts = a->rowners,rank,n,cnt,j; 561 int m,*b_sendj,*garray = a->garray,*lens,*jsendbuf,*a_jsendbuf,*b_jsendbuf; 562 PetscScalar *sendbuf,*recvbuf,*a_sendbuf,*b_sendbuf; 563 564 PetscFunctionBegin; 565 ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr); 566 ierr = MPI_Comm_rank(A->comm,&rank);CHKERRQ(ierr); 567 568 if (scall == MAT_INITIAL_MATRIX) { 569 /* ---------------------------------------------------------------- 570 Tell every processor the number of nonzeros per row 571 */ 572 ierr = PetscMalloc(A->M*sizeof(int),&lens);CHKERRQ(ierr); 573 for (i=a->rstart; i<a->rend; i++) { 574 lens[i] = ad->i[i-a->rstart+1] - ad->i[i-a->rstart] + bd->i[i-a->rstart+1] - bd->i[i-a->rstart]; 575 } 576 sendcount = a->rend - a->rstart; 577 ierr = PetscMalloc(2*size*sizeof(int),&recvcounts);CHKERRQ(ierr); 578 displs = recvcounts + size; 579 for (i=0; i<size; i++) { 580 recvcounts[i] = a->rowners[i+1] - a->rowners[i]; 581 displs[i] = a->rowners[i]; 582 } 583 ierr = MPI_Allgatherv(lens+a->rstart,sendcount,MPI_INT,lens,recvcounts,displs,MPI_INT,A->comm);CHKERRQ(ierr); 584 585 /* --------------------------------------------------------------- 586 Create the sequential matrix of the same type as the local block diagonal 587 */ 588 ierr = MatCreate(PETSC_COMM_SELF,A->M,A->N,PETSC_DETERMINE,PETSC_DETERMINE,&B);CHKERRQ(ierr); 589 ierr = MatSetType(B,a->A->type_name);CHKERRQ(ierr); 590 ierr = MatSeqAIJSetPreallocation(B,0,lens);CHKERRQ(ierr); 591 ierr = PetscMalloc(sizeof(Mat),Bin);CHKERRQ(ierr); 592 **Bin = B; 593 b = (Mat_SeqAIJ *)B->data; 594 595 /*-------------------------------------------------------------------- 596 Copy my part of matrix column indices over 597 */ 598 sendcount = ad->nz + bd->nz; 599 jsendbuf = b->j + b->i[rstarts[rank]]; 600 a_jsendbuf = ad->j; 601 b_jsendbuf = bd->j; 602 n = a->rend - a->rstart; 603 cnt = 0; 604 for (i=0; i<n; i++) { 605 606 /* put in lower diagonal portion */ 607 m = bd->i[i+1] - bd->i[i]; 608 while (m > 0) { 609 /* is it above diagonal (in bd (compressed) numbering) */ 610 if (garray[*b_jsendbuf] > a->rstart + i) break; 611 jsendbuf[cnt++] = garray[*b_jsendbuf++]; 612 m--; 613 } 614 615 /* put in diagonal portion */ 616 for (j=ad->i[i]; j<ad->i[i+1]; j++) { 617 jsendbuf[cnt++] = a->rstart + *a_jsendbuf++; 618 } 619 620 /* put in upper diagonal portion */ 621 while (m-- > 0) { 622 jsendbuf[cnt++] = garray[*b_jsendbuf++]; 623 } 624 } 625 if (cnt != sendcount) SETERRQ2(1,"Corrupted PETSc matrix: nz given %d actual nz %d",sendcount,cnt); 626 627 /*-------------------------------------------------------------------- 628 Gather all column indices to all processors 629 */ 630 for (i=0; i<size; i++) { 631 recvcounts[i] = 0; 632 for (j=a->rowners[i]; j<a->rowners[i+1]; j++) { 633 recvcounts[i] += lens[j]; 634 } 635 } 636 displs[0] = 0; 637 for (i=1; i<size; i++) { 638 displs[i] = displs[i-1] + recvcounts[i-1]; 639 } 640 ierr = MPI_Allgatherv(jsendbuf,sendcount,MPI_INT,b->j,recvcounts,displs,MPI_INT,A->comm);CHKERRQ(ierr); 641 642 /*-------------------------------------------------------------------- 643 Assemble the matrix into useable form (note numerical values not yet set) 644 */ 645 /* set the b->ilen (length of each row) values */ 646 ierr = PetscMemcpy(b->ilen,lens,A->M*sizeof(int));CHKERRQ(ierr); 647 /* set the b->i indices */ 648 b->i[0] = 0; 649 for (i=1; i<=A->M; i++) { 650 b->i[i] = b->i[i-1] + lens[i-1]; 651 } 652 ierr = PetscFree(lens);CHKERRQ(ierr); 653 ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 654 ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 655 656 } else { 657 B = **Bin; 658 b = (Mat_SeqAIJ *)B->data; 659 } 660 661 /*-------------------------------------------------------------------- 662 Copy my part of matrix numerical values into the values location 663 */ 664 sendcount = ad->nz + bd->nz; 665 sendbuf = b->a + b->i[rstarts[rank]]; 666 a_sendbuf = ad->a; 667 b_sendbuf = bd->a; 668 b_sendj = bd->j; 669 n = a->rend - a->rstart; 670 cnt = 0; 671 for (i=0; i<n; i++) { 672 673 /* put in lower diagonal portion */ 674 m = bd->i[i+1] - bd->i[i]; 675 while (m > 0) { 676 /* is it above diagonal (in bd (compressed) numbering) */ 677 if (garray[*b_sendj] > a->rstart + i) break; 678 sendbuf[cnt++] = *b_sendbuf++; 679 m--; 680 b_sendj++; 681 } 682 683 /* put in diagonal portion */ 684 for (j=ad->i[i]; j<ad->i[i+1]; j++) { 685 sendbuf[cnt++] = *a_sendbuf++; 686 } 687 688 /* put in upper diagonal portion */ 689 while (m-- > 0) { 690 sendbuf[cnt++] = *b_sendbuf++; 691 b_sendj++; 692 } 693 } 694 if (cnt != sendcount) SETERRQ2(1,"Corrupted PETSc matrix: nz given %d actual nz %d",sendcount,cnt); 695 696 /* ----------------------------------------------------------------- 697 Gather all numerical values to all processors 698 */ 699 if (!recvcounts) { 700 ierr = PetscMalloc(2*size*sizeof(int),&recvcounts);CHKERRQ(ierr); 701 displs = recvcounts + size; 702 } 703 for (i=0; i<size; i++) { 704 recvcounts[i] = b->i[rstarts[i+1]] - b->i[rstarts[i]]; 705 } 706 displs[0] = 0; 707 for (i=1; i<size; i++) { 708 displs[i] = displs[i-1] + recvcounts[i-1]; 709 } 710 recvbuf = b->a; 711 ierr = MPI_Allgatherv(sendbuf,sendcount,MPIU_SCALAR,recvbuf,recvcounts,displs,MPIU_SCALAR,A->comm);CHKERRQ(ierr); 712 ierr = PetscFree(recvcounts);CHKERRQ(ierr); 713 if (A->symmetric){ 714 ierr = MatSetOption(B,MAT_SYMMETRIC);CHKERRQ(ierr); 715 } else if (A->hermitian) { 716 ierr = MatSetOption(B,MAT_HERMITIAN);CHKERRQ(ierr); 717 } else if (A->structurally_symmetric) { 718 ierr = MatSetOption(B,MAT_STRUCTURALLY_SYMMETRIC);CHKERRQ(ierr); 719 } 720 721 PetscFunctionReturn(0); 722 } 723 724 #undef __FUNCT__ 725 #define __FUNCT__ "MatGetSubMatrices_MPIAIJ" 726 PetscErrorCode MatGetSubMatrices_MPIAIJ(Mat C,int ismax,const IS isrow[],const IS iscol[],MatReuse scall,Mat *submat[]) 727 { 728 PetscErrorCode ierr; 729 int nmax,nstages_local,nstages,i,pos,max_no,nrow,ncol; 730 PetscTruth rowflag,colflag,wantallmatrix = PETSC_FALSE,twantallmatrix; 731 732 PetscFunctionBegin; 733 /* 734 Check for special case each processor gets entire matrix 735 */ 736 if (ismax == 1 && C->M == C->N) { 737 ierr = ISIdentity(*isrow,&rowflag);CHKERRQ(ierr); 738 ierr = ISIdentity(*iscol,&colflag);CHKERRQ(ierr); 739 ierr = ISGetLocalSize(*isrow,&nrow);CHKERRQ(ierr); 740 ierr = ISGetLocalSize(*iscol,&ncol);CHKERRQ(ierr); 741 if (rowflag && colflag && nrow == C->M && ncol == C->N) { 742 wantallmatrix = PETSC_TRUE; 743 ierr = PetscOptionsGetLogical(C->prefix,"-use_fast_submatrix",&wantallmatrix,PETSC_NULL);CHKERRQ(ierr); 744 } 745 } 746 ierr = MPI_Allreduce(&wantallmatrix,&twantallmatrix,1,MPI_INT,MPI_MIN,C->comm);CHKERRQ(ierr); 747 if (twantallmatrix) { 748 ierr = MatGetSubMatrix_MPIAIJ_All(C,scall,submat);CHKERRQ(ierr); 749 PetscFunctionReturn(0); 750 } 751 752 /* Allocate memory to hold all the submatrices */ 753 if (scall != MAT_REUSE_MATRIX) { 754 ierr = PetscMalloc((ismax+1)*sizeof(Mat),submat);CHKERRQ(ierr); 755 } 756 /* Determine the number of stages through which submatrices are done */ 757 nmax = 20*1000000 / (C->N * sizeof(int)); 758 if (!nmax) nmax = 1; 759 nstages_local = ismax/nmax + ((ismax % nmax)?1:0); 760 761 /* Make sure every processor loops through the nstages */ 762 ierr = MPI_Allreduce(&nstages_local,&nstages,1,MPI_INT,MPI_MAX,C->comm);CHKERRQ(ierr); 763 764 for (i=0,pos=0; i<nstages; i++) { 765 if (pos+nmax <= ismax) max_no = nmax; 766 else if (pos == ismax) max_no = 0; 767 else max_no = ismax-pos; 768 ierr = MatGetSubMatrices_MPIAIJ_Local(C,max_no,isrow+pos,iscol+pos,scall,*submat+pos);CHKERRQ(ierr); 769 pos += max_no; 770 } 771 PetscFunctionReturn(0); 772 } 773 /* -------------------------------------------------------------------------*/ 774 #undef __FUNCT__ 775 #define __FUNCT__ "MatGetSubMatrices_MPIAIJ_Local" 776 PetscErrorCode MatGetSubMatrices_MPIAIJ_Local(Mat C,int ismax,const IS isrow[],const IS iscol[],MatReuse scall,Mat *submats) 777 { 778 Mat_MPIAIJ *c = (Mat_MPIAIJ*)C->data; 779 Mat A = c->A; 780 Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*b = (Mat_SeqAIJ*)c->B->data,*mat; 781 int **irow,**icol,*nrow,*ncol,*w1,*w2,*w3,*w4,*rtable,start,end,size; 782 PetscErrorCode ierr; 783 int **sbuf1,**sbuf2,rank,m,i,j,k,l,ct1,ct2,**rbuf1,row,proc; 784 int nrqs,msz,**ptr,idex,*req_size,*ctr,*pa,*tmp,tcol,nrqr; 785 int **rbuf3,*req_source,**sbuf_aj,**rbuf2,max1,max2,**rmap; 786 int **cmap,**lens,is_no,ncols,*cols,mat_i,*mat_j,tmp2,jmax,*irow_i; 787 int len,ctr_j,*sbuf1_j,*sbuf_aj_i,*rbuf1_i,kmax,*cmap_i,*lens_i; 788 int *rmap_i,tag0,tag1,tag2,tag3; 789 MPI_Request *s_waits1,*r_waits1,*s_waits2,*r_waits2,*r_waits3; 790 MPI_Request *r_waits4,*s_waits3,*s_waits4; 791 MPI_Status *r_status1,*r_status2,*s_status1,*s_status3,*s_status2; 792 MPI_Status *r_status3,*r_status4,*s_status4; 793 MPI_Comm comm; 794 PetscScalar **rbuf4,**sbuf_aa,*vals,*mat_a,*sbuf_aa_i; 795 PetscTruth sorted; 796 int *onodes1,*olengths1; 797 798 PetscFunctionBegin; 799 comm = C->comm; 800 tag0 = C->tag; 801 size = c->size; 802 rank = c->rank; 803 m = C->M; 804 805 /* Get some new tags to keep the communication clean */ 806 ierr = PetscObjectGetNewTag((PetscObject)C,&tag1);CHKERRQ(ierr); 807 ierr = PetscObjectGetNewTag((PetscObject)C,&tag2);CHKERRQ(ierr); 808 ierr = PetscObjectGetNewTag((PetscObject)C,&tag3);CHKERRQ(ierr); 809 810 /* Check if the col indices are sorted */ 811 for (i=0; i<ismax; i++) { 812 ierr = ISSorted(isrow[i],&sorted);CHKERRQ(ierr); 813 if (!sorted) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"ISrow is not sorted"); 814 ierr = ISSorted(iscol[i],&sorted);CHKERRQ(ierr); 815 /* if (!sorted) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"IScol is not sorted"); */ 816 } 817 818 len = (2*ismax+1)*(sizeof(int*)+ sizeof(int)) + (m+1)*sizeof(int); 819 ierr = PetscMalloc(len,&irow);CHKERRQ(ierr); 820 icol = irow + ismax; 821 nrow = (int*)(icol + ismax); 822 ncol = nrow + ismax; 823 rtable = ncol + ismax; 824 825 for (i=0; i<ismax; i++) { 826 ierr = ISGetIndices(isrow[i],&irow[i]);CHKERRQ(ierr); 827 ierr = ISGetIndices(iscol[i],&icol[i]);CHKERRQ(ierr); 828 ierr = ISGetLocalSize(isrow[i],&nrow[i]);CHKERRQ(ierr); 829 ierr = ISGetLocalSize(iscol[i],&ncol[i]);CHKERRQ(ierr); 830 } 831 832 /* Create hash table for the mapping :row -> proc*/ 833 for (i=0,j=0; i<size; i++) { 834 jmax = c->rowners[i+1]; 835 for (; j<jmax; j++) { 836 rtable[j] = i; 837 } 838 } 839 840 /* evaluate communication - mesg to who, length of mesg, and buffer space 841 required. Based on this, buffers are allocated, and data copied into them*/ 842 ierr = PetscMalloc(size*4*sizeof(int),&w1);CHKERRQ(ierr); /* mesg size */ 843 w2 = w1 + size; /* if w2[i] marked, then a message to proc i*/ 844 w3 = w2 + size; /* no of IS that needs to be sent to proc i */ 845 w4 = w3 + size; /* temp work space used in determining w1, w2, w3 */ 846 ierr = PetscMemzero(w1,size*3*sizeof(int));CHKERRQ(ierr); /* initialize work vector*/ 847 for (i=0; i<ismax; i++) { 848 ierr = PetscMemzero(w4,size*sizeof(int));CHKERRQ(ierr); /* initialize work vector*/ 849 jmax = nrow[i]; 850 irow_i = irow[i]; 851 for (j=0; j<jmax; j++) { 852 row = irow_i[j]; 853 proc = rtable[row]; 854 w4[proc]++; 855 } 856 for (j=0; j<size; j++) { 857 if (w4[j]) { w1[j] += w4[j]; w3[j]++;} 858 } 859 } 860 861 nrqs = 0; /* no of outgoing messages */ 862 msz = 0; /* total mesg length (for all procs) */ 863 w1[rank] = 0; /* no mesg sent to self */ 864 w3[rank] = 0; 865 for (i=0; i<size; i++) { 866 if (w1[i]) { w2[i] = 1; nrqs++;} /* there exists a message to proc i */ 867 } 868 ierr = PetscMalloc((nrqs+1)*sizeof(int),&pa);CHKERRQ(ierr); /*(proc -array)*/ 869 for (i=0,j=0; i<size; i++) { 870 if (w1[i]) { pa[j] = i; j++; } 871 } 872 873 /* Each message would have a header = 1 + 2*(no of IS) + data */ 874 for (i=0; i<nrqs; i++) { 875 j = pa[i]; 876 w1[j] += w2[j] + 2* w3[j]; 877 msz += w1[j]; 878 } 879 880 /* Determine the number of messages to expect, their lengths, from from-ids */ 881 ierr = PetscGatherNumberOfMessages(comm,w2,w1,&nrqr);CHKERRQ(ierr); 882 ierr = PetscGatherMessageLengths(comm,nrqs,nrqr,w1,&onodes1,&olengths1);CHKERRQ(ierr); 883 884 /* Now post the Irecvs corresponding to these messages */ 885 ierr = PetscPostIrecvInt(comm,tag0,nrqr,onodes1,olengths1,&rbuf1,&r_waits1);CHKERRQ(ierr); 886 887 ierr = PetscFree(onodes1);CHKERRQ(ierr); 888 ierr = PetscFree(olengths1);CHKERRQ(ierr); 889 890 /* Allocate Memory for outgoing messages */ 891 len = 2*size*sizeof(int*) + 2*msz*sizeof(int) + size*sizeof(int); 892 ierr = PetscMalloc(len,&sbuf1);CHKERRQ(ierr); 893 ptr = sbuf1 + size; /* Pointers to the data in outgoing buffers */ 894 ierr = PetscMemzero(sbuf1,2*size*sizeof(int*));CHKERRQ(ierr); 895 /* allocate memory for outgoing data + buf to receive the first reply */ 896 tmp = (int*)(ptr + size); 897 ctr = tmp + 2*msz; 898 899 { 900 int *iptr = tmp,ict = 0; 901 for (i=0; i<nrqs; i++) { 902 j = pa[i]; 903 iptr += ict; 904 sbuf1[j] = iptr; 905 ict = w1[j]; 906 } 907 } 908 909 /* Form the outgoing messages */ 910 /* Initialize the header space */ 911 for (i=0; i<nrqs; i++) { 912 j = pa[i]; 913 sbuf1[j][0] = 0; 914 ierr = PetscMemzero(sbuf1[j]+1,2*w3[j]*sizeof(int));CHKERRQ(ierr); 915 ptr[j] = sbuf1[j] + 2*w3[j] + 1; 916 } 917 918 /* Parse the isrow and copy data into outbuf */ 919 for (i=0; i<ismax; i++) { 920 ierr = PetscMemzero(ctr,size*sizeof(int));CHKERRQ(ierr); 921 irow_i = irow[i]; 922 jmax = nrow[i]; 923 for (j=0; j<jmax; j++) { /* parse the indices of each IS */ 924 row = irow_i[j]; 925 proc = rtable[row]; 926 if (proc != rank) { /* copy to the outgoing buf*/ 927 ctr[proc]++; 928 *ptr[proc] = row; 929 ptr[proc]++; 930 } 931 } 932 /* Update the headers for the current IS */ 933 for (j=0; j<size; j++) { /* Can Optimise this loop too */ 934 if ((ctr_j = ctr[j])) { 935 sbuf1_j = sbuf1[j]; 936 k = ++sbuf1_j[0]; 937 sbuf1_j[2*k] = ctr_j; 938 sbuf1_j[2*k-1] = i; 939 } 940 } 941 } 942 943 /* Now post the sends */ 944 ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Request),&s_waits1);CHKERRQ(ierr); 945 for (i=0; i<nrqs; ++i) { 946 j = pa[i]; 947 ierr = MPI_Isend(sbuf1[j],w1[j],MPI_INT,j,tag0,comm,s_waits1+i);CHKERRQ(ierr); 948 } 949 950 /* Post Receives to capture the buffer size */ 951 ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Request),&r_waits2);CHKERRQ(ierr); 952 ierr = PetscMalloc((nrqs+1)*sizeof(int*),&rbuf2);CHKERRQ(ierr); 953 rbuf2[0] = tmp + msz; 954 for (i=1; i<nrqs; ++i) { 955 rbuf2[i] = rbuf2[i-1]+w1[pa[i-1]]; 956 } 957 for (i=0; i<nrqs; ++i) { 958 j = pa[i]; 959 ierr = MPI_Irecv(rbuf2[i],w1[j],MPI_INT,j,tag1,comm,r_waits2+i);CHKERRQ(ierr); 960 } 961 962 /* Send to other procs the buf size they should allocate */ 963 964 965 /* Receive messages*/ 966 ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Request),&s_waits2);CHKERRQ(ierr); 967 ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Status),&r_status1);CHKERRQ(ierr); 968 len = 2*nrqr*sizeof(int) + (nrqr+1)*sizeof(int*); 969 ierr = PetscMalloc(len,&sbuf2);CHKERRQ(ierr); 970 req_size = (int*)(sbuf2 + nrqr); 971 req_source = req_size + nrqr; 972 973 { 974 Mat_SeqAIJ *sA = (Mat_SeqAIJ*)c->A->data,*sB = (Mat_SeqAIJ*)c->B->data; 975 int *sAi = sA->i,*sBi = sB->i,id,rstart = c->rstart; 976 int *sbuf2_i; 977 978 for (i=0; i<nrqr; ++i) { 979 ierr = MPI_Waitany(nrqr,r_waits1,&idex,r_status1+i);CHKERRQ(ierr); 980 req_size[idex] = 0; 981 rbuf1_i = rbuf1[idex]; 982 start = 2*rbuf1_i[0] + 1; 983 ierr = MPI_Get_count(r_status1+i,MPI_INT,&end);CHKERRQ(ierr); 984 ierr = PetscMalloc((end+1)*sizeof(int),&sbuf2[idex]);CHKERRQ(ierr); 985 sbuf2_i = sbuf2[idex]; 986 for (j=start; j<end; j++) { 987 id = rbuf1_i[j] - rstart; 988 ncols = sAi[id+1] - sAi[id] + sBi[id+1] - sBi[id]; 989 sbuf2_i[j] = ncols; 990 req_size[idex] += ncols; 991 } 992 req_source[idex] = r_status1[i].MPI_SOURCE; 993 /* form the header */ 994 sbuf2_i[0] = req_size[idex]; 995 for (j=1; j<start; j++) { sbuf2_i[j] = rbuf1_i[j]; } 996 ierr = MPI_Isend(sbuf2_i,end,MPI_INT,req_source[idex],tag1,comm,s_waits2+i);CHKERRQ(ierr); 997 } 998 } 999 ierr = PetscFree(r_status1);CHKERRQ(ierr); 1000 ierr = PetscFree(r_waits1);CHKERRQ(ierr); 1001 1002 /* recv buffer sizes */ 1003 /* Receive messages*/ 1004 1005 ierr = PetscMalloc((nrqs+1)*sizeof(int*),&rbuf3);CHKERRQ(ierr); 1006 ierr = PetscMalloc((nrqs+1)*sizeof(PetscScalar*),&rbuf4);CHKERRQ(ierr); 1007 ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Request),&r_waits3);CHKERRQ(ierr); 1008 ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Request),&r_waits4);CHKERRQ(ierr); 1009 ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Status),&r_status2);CHKERRQ(ierr); 1010 1011 for (i=0; i<nrqs; ++i) { 1012 ierr = MPI_Waitany(nrqs,r_waits2,&idex,r_status2+i);CHKERRQ(ierr); 1013 ierr = PetscMalloc((rbuf2[idex][0]+1)*sizeof(int),&rbuf3[idex]);CHKERRQ(ierr); 1014 ierr = PetscMalloc((rbuf2[idex][0]+1)*sizeof(PetscScalar),&rbuf4[idex]);CHKERRQ(ierr); 1015 ierr = MPI_Irecv(rbuf3[idex],rbuf2[idex][0],MPI_INT,r_status2[i].MPI_SOURCE,tag2,comm,r_waits3+idex);CHKERRQ(ierr); 1016 ierr = MPI_Irecv(rbuf4[idex],rbuf2[idex][0],MPIU_SCALAR,r_status2[i].MPI_SOURCE,tag3,comm,r_waits4+idex);CHKERRQ(ierr); 1017 } 1018 ierr = PetscFree(r_status2);CHKERRQ(ierr); 1019 ierr = PetscFree(r_waits2);CHKERRQ(ierr); 1020 1021 /* Wait on sends1 and sends2 */ 1022 ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Status),&s_status1);CHKERRQ(ierr); 1023 ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Status),&s_status2);CHKERRQ(ierr); 1024 1025 ierr = MPI_Waitall(nrqs,s_waits1,s_status1);CHKERRQ(ierr); 1026 ierr = MPI_Waitall(nrqr,s_waits2,s_status2);CHKERRQ(ierr); 1027 ierr = PetscFree(s_status1);CHKERRQ(ierr); 1028 ierr = PetscFree(s_status2);CHKERRQ(ierr); 1029 ierr = PetscFree(s_waits1);CHKERRQ(ierr); 1030 ierr = PetscFree(s_waits2);CHKERRQ(ierr); 1031 1032 /* Now allocate buffers for a->j, and send them off */ 1033 ierr = PetscMalloc((nrqr+1)*sizeof(int*),&sbuf_aj);CHKERRQ(ierr); 1034 for (i=0,j=0; i<nrqr; i++) j += req_size[i]; 1035 ierr = PetscMalloc((j+1)*sizeof(int),&sbuf_aj[0]);CHKERRQ(ierr); 1036 for (i=1; i<nrqr; i++) sbuf_aj[i] = sbuf_aj[i-1] + req_size[i-1]; 1037 1038 ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Request),&s_waits3);CHKERRQ(ierr); 1039 { 1040 int nzA,nzB,*a_i = a->i,*b_i = b->i,imark; 1041 int *cworkA,*cworkB,cstart = c->cstart,rstart = c->rstart,*bmap = c->garray; 1042 int *a_j = a->j,*b_j = b->j,ctmp; 1043 1044 for (i=0; i<nrqr; i++) { 1045 rbuf1_i = rbuf1[i]; 1046 sbuf_aj_i = sbuf_aj[i]; 1047 ct1 = 2*rbuf1_i[0] + 1; 1048 ct2 = 0; 1049 for (j=1,max1=rbuf1_i[0]; j<=max1; j++) { 1050 kmax = rbuf1[i][2*j]; 1051 for (k=0; k<kmax; k++,ct1++) { 1052 row = rbuf1_i[ct1] - rstart; 1053 nzA = a_i[row+1] - a_i[row]; nzB = b_i[row+1] - b_i[row]; 1054 ncols = nzA + nzB; 1055 cworkA = a_j + a_i[row]; cworkB = b_j + b_i[row]; 1056 1057 /* load the column indices for this row into cols*/ 1058 cols = sbuf_aj_i + ct2; 1059 1060 for (l=0; l<nzB; l++) { 1061 if ((ctmp = bmap[cworkB[l]]) < cstart) cols[l] = ctmp; 1062 else break; 1063 } 1064 imark = l; 1065 for (l=0; l<nzA; l++) cols[imark+l] = cstart + cworkA[l]; 1066 for (l=imark; l<nzB; l++) cols[nzA+l] = bmap[cworkB[l]]; 1067 1068 ct2 += ncols; 1069 } 1070 } 1071 ierr = MPI_Isend(sbuf_aj_i,req_size[i],MPI_INT,req_source[i],tag2,comm,s_waits3+i);CHKERRQ(ierr); 1072 } 1073 } 1074 ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Status),&r_status3);CHKERRQ(ierr); 1075 ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Status),&s_status3);CHKERRQ(ierr); 1076 1077 /* Allocate buffers for a->a, and send them off */ 1078 ierr = PetscMalloc((nrqr+1)*sizeof(PetscScalar*),&sbuf_aa);CHKERRQ(ierr); 1079 for (i=0,j=0; i<nrqr; i++) j += req_size[i]; 1080 ierr = PetscMalloc((j+1)*sizeof(PetscScalar),&sbuf_aa[0]);CHKERRQ(ierr); 1081 for (i=1; i<nrqr; i++) sbuf_aa[i] = sbuf_aa[i-1] + req_size[i-1]; 1082 1083 ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Request),&s_waits4);CHKERRQ(ierr); 1084 { 1085 int nzA,nzB,*a_i = a->i,*b_i = b->i, *cworkB,imark; 1086 int cstart = c->cstart,rstart = c->rstart,*bmap = c->garray; 1087 int *b_j = b->j; 1088 PetscScalar *vworkA,*vworkB,*a_a = a->a,*b_a = b->a; 1089 1090 for (i=0; i<nrqr; i++) { 1091 rbuf1_i = rbuf1[i]; 1092 sbuf_aa_i = sbuf_aa[i]; 1093 ct1 = 2*rbuf1_i[0]+1; 1094 ct2 = 0; 1095 for (j=1,max1=rbuf1_i[0]; j<=max1; j++) { 1096 kmax = rbuf1_i[2*j]; 1097 for (k=0; k<kmax; k++,ct1++) { 1098 row = rbuf1_i[ct1] - rstart; 1099 nzA = a_i[row+1] - a_i[row]; nzB = b_i[row+1] - b_i[row]; 1100 ncols = nzA + nzB; 1101 cworkB = b_j + b_i[row]; 1102 vworkA = a_a + a_i[row]; 1103 vworkB = b_a + b_i[row]; 1104 1105 /* load the column values for this row into vals*/ 1106 vals = sbuf_aa_i+ct2; 1107 1108 for (l=0; l<nzB; l++) { 1109 if ((bmap[cworkB[l]]) < cstart) vals[l] = vworkB[l]; 1110 else break; 1111 } 1112 imark = l; 1113 for (l=0; l<nzA; l++) vals[imark+l] = vworkA[l]; 1114 for (l=imark; l<nzB; l++) vals[nzA+l] = vworkB[l]; 1115 1116 ct2 += ncols; 1117 } 1118 } 1119 ierr = MPI_Isend(sbuf_aa_i,req_size[i],MPIU_SCALAR,req_source[i],tag3,comm,s_waits4+i);CHKERRQ(ierr); 1120 } 1121 } 1122 ierr = PetscMalloc((nrqs+1)*sizeof(MPI_Status),&r_status4);CHKERRQ(ierr); 1123 ierr = PetscMalloc((nrqr+1)*sizeof(MPI_Status),&s_status4);CHKERRQ(ierr); 1124 ierr = PetscFree(rbuf1);CHKERRQ(ierr); 1125 1126 /* Form the matrix */ 1127 /* create col map */ 1128 { 1129 int *icol_i; 1130 1131 len = (1+ismax)*sizeof(int*)+ ismax*C->N*sizeof(int); 1132 ierr = PetscMalloc(len,&cmap);CHKERRQ(ierr); 1133 cmap[0] = (int*)(cmap + ismax); 1134 ierr = PetscMemzero(cmap[0],(1+ismax*C->N)*sizeof(int));CHKERRQ(ierr); 1135 for (i=1; i<ismax; i++) { cmap[i] = cmap[i-1] + C->N; } 1136 for (i=0; i<ismax; i++) { 1137 jmax = ncol[i]; 1138 icol_i = icol[i]; 1139 cmap_i = cmap[i]; 1140 for (j=0; j<jmax; j++) { 1141 cmap_i[icol_i[j]] = j+1; 1142 } 1143 } 1144 } 1145 1146 /* Create lens which is required for MatCreate... */ 1147 for (i=0,j=0; i<ismax; i++) { j += nrow[i]; } 1148 len = (1+ismax)*sizeof(int*)+ j*sizeof(int); 1149 ierr = PetscMalloc(len,&lens);CHKERRQ(ierr); 1150 lens[0] = (int*)(lens + ismax); 1151 ierr = PetscMemzero(lens[0],j*sizeof(int));CHKERRQ(ierr); 1152 for (i=1; i<ismax; i++) { lens[i] = lens[i-1] + nrow[i-1]; } 1153 1154 /* Update lens from local data */ 1155 for (i=0; i<ismax; i++) { 1156 jmax = nrow[i]; 1157 cmap_i = cmap[i]; 1158 irow_i = irow[i]; 1159 lens_i = lens[i]; 1160 for (j=0; j<jmax; j++) { 1161 row = irow_i[j]; 1162 proc = rtable[row]; 1163 if (proc == rank) { 1164 ierr = MatGetRow_MPIAIJ(C,row,&ncols,&cols,0);CHKERRQ(ierr); 1165 for (k=0; k<ncols; k++) { 1166 if (cmap_i[cols[k]]) { lens_i[j]++;} 1167 } 1168 ierr = MatRestoreRow_MPIAIJ(C,row,&ncols,&cols,0);CHKERRQ(ierr); 1169 } 1170 } 1171 } 1172 1173 /* Create row map*/ 1174 len = (1+ismax)*sizeof(int*)+ ismax*C->M*sizeof(int); 1175 ierr = PetscMalloc(len,&rmap);CHKERRQ(ierr); 1176 rmap[0] = (int*)(rmap + ismax); 1177 ierr = PetscMemzero(rmap[0],ismax*C->M*sizeof(int));CHKERRQ(ierr); 1178 for (i=1; i<ismax; i++) { rmap[i] = rmap[i-1] + C->M;} 1179 for (i=0; i<ismax; i++) { 1180 rmap_i = rmap[i]; 1181 irow_i = irow[i]; 1182 jmax = nrow[i]; 1183 for (j=0; j<jmax; j++) { 1184 rmap_i[irow_i[j]] = j; 1185 } 1186 } 1187 1188 /* Update lens from offproc data */ 1189 { 1190 int *rbuf2_i,*rbuf3_i,*sbuf1_i; 1191 1192 for (tmp2=0; tmp2<nrqs; tmp2++) { 1193 ierr = MPI_Waitany(nrqs,r_waits3,&i,r_status3+tmp2);CHKERRQ(ierr); 1194 idex = pa[i]; 1195 sbuf1_i = sbuf1[idex]; 1196 jmax = sbuf1_i[0]; 1197 ct1 = 2*jmax+1; 1198 ct2 = 0; 1199 rbuf2_i = rbuf2[i]; 1200 rbuf3_i = rbuf3[i]; 1201 for (j=1; j<=jmax; j++) { 1202 is_no = sbuf1_i[2*j-1]; 1203 max1 = sbuf1_i[2*j]; 1204 lens_i = lens[is_no]; 1205 cmap_i = cmap[is_no]; 1206 rmap_i = rmap[is_no]; 1207 for (k=0; k<max1; k++,ct1++) { 1208 row = rmap_i[sbuf1_i[ct1]]; /* the val in the new matrix to be */ 1209 max2 = rbuf2_i[ct1]; 1210 for (l=0; l<max2; l++,ct2++) { 1211 if (cmap_i[rbuf3_i[ct2]]) { 1212 lens_i[row]++; 1213 } 1214 } 1215 } 1216 } 1217 } 1218 } 1219 ierr = PetscFree(r_status3);CHKERRQ(ierr); 1220 ierr = PetscFree(r_waits3);CHKERRQ(ierr); 1221 ierr = MPI_Waitall(nrqr,s_waits3,s_status3);CHKERRQ(ierr); 1222 ierr = PetscFree(s_status3);CHKERRQ(ierr); 1223 ierr = PetscFree(s_waits3);CHKERRQ(ierr); 1224 1225 /* Create the submatrices */ 1226 if (scall == MAT_REUSE_MATRIX) { 1227 PetscTruth flag; 1228 1229 /* 1230 Assumes new rows are same length as the old rows,hence bug! 1231 */ 1232 for (i=0; i<ismax; i++) { 1233 mat = (Mat_SeqAIJ *)(submats[i]->data); 1234 if ((submats[i]->m != nrow[i]) || (submats[i]->n != ncol[i])) { 1235 SETERRQ(PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong size"); 1236 } 1237 ierr = PetscMemcmp(mat->ilen,lens[i],submats[i]->m*sizeof(int),&flag);CHKERRQ(ierr); 1238 if (flag == PETSC_FALSE) { 1239 SETERRQ(PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong no of nonzeros"); 1240 } 1241 /* Initial matrix as if empty */ 1242 ierr = PetscMemzero(mat->ilen,submats[i]->m*sizeof(int));CHKERRQ(ierr); 1243 submats[i]->factor = C->factor; 1244 } 1245 } else { 1246 for (i=0; i<ismax; i++) { 1247 ierr = MatCreate(PETSC_COMM_SELF,nrow[i],ncol[i],PETSC_DETERMINE,PETSC_DETERMINE,submats+i);CHKERRQ(ierr); 1248 ierr = MatSetType(submats[i],A->type_name);CHKERRQ(ierr); 1249 ierr = MatSeqAIJSetPreallocation(submats[i],0,lens[i]);CHKERRQ(ierr); 1250 } 1251 } 1252 1253 /* Assemble the matrices */ 1254 /* First assemble the local rows */ 1255 { 1256 int ilen_row,*imat_ilen,*imat_j,*imat_i,old_row; 1257 PetscScalar *imat_a; 1258 1259 for (i=0; i<ismax; i++) { 1260 mat = (Mat_SeqAIJ*)submats[i]->data; 1261 imat_ilen = mat->ilen; 1262 imat_j = mat->j; 1263 imat_i = mat->i; 1264 imat_a = mat->a; 1265 cmap_i = cmap[i]; 1266 rmap_i = rmap[i]; 1267 irow_i = irow[i]; 1268 jmax = nrow[i]; 1269 for (j=0; j<jmax; j++) { 1270 row = irow_i[j]; 1271 proc = rtable[row]; 1272 if (proc == rank) { 1273 old_row = row; 1274 row = rmap_i[row]; 1275 ilen_row = imat_ilen[row]; 1276 ierr = MatGetRow_MPIAIJ(C,old_row,&ncols,&cols,&vals);CHKERRQ(ierr); 1277 mat_i = imat_i[row] ; 1278 mat_a = imat_a + mat_i; 1279 mat_j = imat_j + mat_i; 1280 for (k=0; k<ncols; k++) { 1281 if ((tcol = cmap_i[cols[k]])) { 1282 *mat_j++ = tcol - 1; 1283 *mat_a++ = vals[k]; 1284 ilen_row++; 1285 } 1286 } 1287 ierr = MatRestoreRow_MPIAIJ(C,old_row,&ncols,&cols,&vals);CHKERRQ(ierr); 1288 imat_ilen[row] = ilen_row; 1289 } 1290 } 1291 } 1292 } 1293 1294 /* Now assemble the off proc rows*/ 1295 { 1296 int *sbuf1_i,*rbuf2_i,*rbuf3_i,*imat_ilen,ilen; 1297 int *imat_j,*imat_i; 1298 PetscScalar *imat_a,*rbuf4_i; 1299 1300 for (tmp2=0; tmp2<nrqs; tmp2++) { 1301 ierr = MPI_Waitany(nrqs,r_waits4,&i,r_status4+tmp2);CHKERRQ(ierr); 1302 idex = pa[i]; 1303 sbuf1_i = sbuf1[idex]; 1304 jmax = sbuf1_i[0]; 1305 ct1 = 2*jmax + 1; 1306 ct2 = 0; 1307 rbuf2_i = rbuf2[i]; 1308 rbuf3_i = rbuf3[i]; 1309 rbuf4_i = rbuf4[i]; 1310 for (j=1; j<=jmax; j++) { 1311 is_no = sbuf1_i[2*j-1]; 1312 rmap_i = rmap[is_no]; 1313 cmap_i = cmap[is_no]; 1314 mat = (Mat_SeqAIJ*)submats[is_no]->data; 1315 imat_ilen = mat->ilen; 1316 imat_j = mat->j; 1317 imat_i = mat->i; 1318 imat_a = mat->a; 1319 max1 = sbuf1_i[2*j]; 1320 for (k=0; k<max1; k++,ct1++) { 1321 row = sbuf1_i[ct1]; 1322 row = rmap_i[row]; 1323 ilen = imat_ilen[row]; 1324 mat_i = imat_i[row] ; 1325 mat_a = imat_a + mat_i; 1326 mat_j = imat_j + mat_i; 1327 max2 = rbuf2_i[ct1]; 1328 for (l=0; l<max2; l++,ct2++) { 1329 if ((tcol = cmap_i[rbuf3_i[ct2]])) { 1330 *mat_j++ = tcol - 1; 1331 *mat_a++ = rbuf4_i[ct2]; 1332 ilen++; 1333 } 1334 } 1335 imat_ilen[row] = ilen; 1336 } 1337 } 1338 } 1339 } 1340 ierr = PetscFree(r_status4);CHKERRQ(ierr); 1341 ierr = PetscFree(r_waits4);CHKERRQ(ierr); 1342 ierr = MPI_Waitall(nrqr,s_waits4,s_status4);CHKERRQ(ierr); 1343 ierr = PetscFree(s_waits4);CHKERRQ(ierr); 1344 ierr = PetscFree(s_status4);CHKERRQ(ierr); 1345 1346 /* Restore the indices */ 1347 for (i=0; i<ismax; i++) { 1348 ierr = ISRestoreIndices(isrow[i],irow+i);CHKERRQ(ierr); 1349 ierr = ISRestoreIndices(iscol[i],icol+i);CHKERRQ(ierr); 1350 } 1351 1352 /* Destroy allocated memory */ 1353 ierr = PetscFree(irow);CHKERRQ(ierr); 1354 ierr = PetscFree(w1);CHKERRQ(ierr); 1355 ierr = PetscFree(pa);CHKERRQ(ierr); 1356 1357 ierr = PetscFree(sbuf1);CHKERRQ(ierr); 1358 ierr = PetscFree(rbuf2);CHKERRQ(ierr); 1359 for (i=0; i<nrqr; ++i) { 1360 ierr = PetscFree(sbuf2[i]);CHKERRQ(ierr); 1361 } 1362 for (i=0; i<nrqs; ++i) { 1363 ierr = PetscFree(rbuf3[i]);CHKERRQ(ierr); 1364 ierr = PetscFree(rbuf4[i]);CHKERRQ(ierr); 1365 } 1366 1367 ierr = PetscFree(sbuf2);CHKERRQ(ierr); 1368 ierr = PetscFree(rbuf3);CHKERRQ(ierr); 1369 ierr = PetscFree(rbuf4);CHKERRQ(ierr); 1370 ierr = PetscFree(sbuf_aj[0]);CHKERRQ(ierr); 1371 ierr = PetscFree(sbuf_aj);CHKERRQ(ierr); 1372 ierr = PetscFree(sbuf_aa[0]);CHKERRQ(ierr); 1373 ierr = PetscFree(sbuf_aa);CHKERRQ(ierr); 1374 1375 ierr = PetscFree(cmap);CHKERRQ(ierr); 1376 ierr = PetscFree(rmap);CHKERRQ(ierr); 1377 ierr = PetscFree(lens);CHKERRQ(ierr); 1378 1379 for (i=0; i<ismax; i++) { 1380 ierr = MatAssemblyBegin(submats[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1381 ierr = MatAssemblyEnd(submats[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1382 } 1383 PetscFunctionReturn(0); 1384 } 1385 1386