1 2 /* 3 Routines to compute overlapping regions of a parallel MPI matrix 4 and to find submatrices that were shared across processors. 5 */ 6 #include <../src/mat/impls/baij/mpi/mpibaij.h> 7 #include <petscbt.h> 8 9 static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Local(Mat,PetscInt,char**,PetscInt*,PetscInt**); 10 static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Receive(Mat,PetscInt,PetscInt**,PetscInt**,PetscInt*); 11 extern PetscErrorCode MatGetRow_MPIBAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**); 12 extern PetscErrorCode MatRestoreRow_MPIBAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**); 13 14 PetscErrorCode MatIncreaseOverlap_MPIBAIJ(Mat C,PetscInt imax,IS is[],PetscInt ov) 15 { 16 PetscErrorCode ierr; 17 PetscInt i,N=C->cmap->N, bs=C->rmap->bs; 18 IS *is_new; 19 20 PetscFunctionBegin; 21 ierr = PetscMalloc1(imax,&is_new);CHKERRQ(ierr); 22 /* Convert the indices into block format */ 23 ierr = ISCompressIndicesGeneral(N,C->rmap->n,bs,imax,is,is_new);CHKERRQ(ierr); 24 if (ov < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap specified\n"); 25 for (i=0; i<ov; ++i) { 26 ierr = MatIncreaseOverlap_MPIBAIJ_Once(C,imax,is_new);CHKERRQ(ierr); 27 } 28 for (i=0; i<imax; i++) {ierr = ISDestroy(&is[i]);CHKERRQ(ierr);} 29 ierr = ISExpandIndicesGeneral(N,N,bs,imax,is_new,is);CHKERRQ(ierr); 30 for (i=0; i<imax; i++) {ierr = ISDestroy(&is_new[i]);CHKERRQ(ierr);} 31 ierr = PetscFree(is_new);CHKERRQ(ierr); 32 PetscFunctionReturn(0); 33 } 34 35 /* 36 Sample message format: 37 If a processor A wants processor B to process some elements corresponding 38 to index sets is[1], is[5] 39 mesg [0] = 2 (no of index sets in the mesg) 40 ----------- 41 mesg [1] = 1 => is[1] 42 mesg [2] = sizeof(is[1]); 43 ----------- 44 mesg [5] = 5 => is[5] 45 mesg [6] = sizeof(is[5]); 46 ----------- 47 mesg [7] 48 mesg [n] data(is[1]) 49 ----------- 50 mesg[n+1] 51 mesg[m] data(is[5]) 52 ----------- 53 54 Notes: 55 nrqs - no of requests sent (or to be sent out) 56 nrqr - no of requests recieved (which have to be or which have been processed 57 */ 58 PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Once(Mat C,PetscInt imax,IS is[]) 59 { 60 Mat_MPIBAIJ *c = (Mat_MPIBAIJ*)C->data; 61 const PetscInt **idx,*idx_i; 62 PetscInt *n,*w3,*w4,**data,len; 63 PetscErrorCode ierr; 64 PetscMPIInt size,rank,tag1,tag2,*w2,*w1,nrqr; 65 PetscInt Mbs,i,j,k,**rbuf,row,proc=-1,nrqs,msz,**outdat,**ptr; 66 PetscInt *ctr,*pa,*tmp,*isz,*isz1,**xdata,**rbuf2,*d_p; 67 PetscMPIInt *onodes1,*olengths1,*onodes2,*olengths2; 68 PetscBT *table; 69 MPI_Comm comm; 70 MPI_Request *s_waits1,*r_waits1,*s_waits2,*r_waits2; 71 MPI_Status *s_status,*recv_status; 72 char *t_p; 73 74 PetscFunctionBegin; 75 ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr); 76 size = c->size; 77 rank = c->rank; 78 Mbs = c->Mbs; 79 80 ierr = PetscObjectGetNewTag((PetscObject)C,&tag1);CHKERRQ(ierr); 81 ierr = PetscObjectGetNewTag((PetscObject)C,&tag2);CHKERRQ(ierr); 82 83 ierr = PetscMalloc2(imax+1,&idx,imax,&n);CHKERRQ(ierr); 84 85 for (i=0; i<imax; i++) { 86 ierr = ISGetIndices(is[i],&idx[i]);CHKERRQ(ierr); 87 ierr = ISGetLocalSize(is[i],&n[i]);CHKERRQ(ierr); 88 } 89 90 /* evaluate communication - mesg to who,length of mesg, and buffer space 91 required. Based on this, buffers are allocated, and data copied into them*/ 92 ierr = PetscCalloc4(size,&w1,size,&w2,size,&w3,size,&w4);CHKERRQ(ierr); 93 for (i=0; i<imax; i++) { 94 ierr = PetscMemzero(w4,size*sizeof(PetscInt));CHKERRQ(ierr); /* initialise work vector*/ 95 idx_i = idx[i]; 96 len = n[i]; 97 for (j=0; j<len; j++) { 98 row = idx_i[j]; 99 if (row < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Index set cannot have negative entries"); 100 ierr = PetscLayoutFindOwner(C->rmap,row*C->rmap->bs,&proc);CHKERRQ(ierr); 101 w4[proc]++; 102 } 103 for (j=0; j<size; j++) { 104 if (w4[j]) { w1[j] += w4[j]; w3[j]++;} 105 } 106 } 107 108 nrqs = 0; /* no of outgoing messages */ 109 msz = 0; /* total mesg length (for all proc */ 110 w1[rank] = 0; /* no mesg sent to itself */ 111 w3[rank] = 0; 112 for (i=0; i<size; i++) { 113 if (w1[i]) {w2[i] = 1; nrqs++;} /* there exists a message to proc i */ 114 } 115 /* pa - is list of processors to communicate with */ 116 ierr = PetscMalloc1(nrqs+1,&pa);CHKERRQ(ierr); 117 for (i=0,j=0; i<size; i++) { 118 if (w1[i]) {pa[j] = i; j++;} 119 } 120 121 /* Each message would have a header = 1 + 2*(no of IS) + data */ 122 for (i=0; i<nrqs; i++) { 123 j = pa[i]; 124 w1[j] += w2[j] + 2*w3[j]; 125 msz += w1[j]; 126 } 127 128 /* Determine the number of messages to expect, their lengths, from from-ids */ 129 ierr = PetscGatherNumberOfMessages(comm,w2,w1,&nrqr);CHKERRQ(ierr); 130 ierr = PetscGatherMessageLengths(comm,nrqs,nrqr,w1,&onodes1,&olengths1);CHKERRQ(ierr); 131 132 /* Now post the Irecvs corresponding to these messages */ 133 ierr = PetscPostIrecvInt(comm,tag1,nrqr,onodes1,olengths1,&rbuf,&r_waits1);CHKERRQ(ierr); 134 135 /* Allocate Memory for outgoing messages */ 136 ierr = PetscMalloc4(size,&outdat,size,&ptr,msz,&tmp,size,&ctr);CHKERRQ(ierr); 137 ierr = PetscMemzero(outdat,size*sizeof(PetscInt*));CHKERRQ(ierr); 138 ierr = PetscMemzero(ptr,size*sizeof(PetscInt*));CHKERRQ(ierr); 139 { 140 PetscInt *iptr = tmp,ict = 0; 141 for (i=0; i<nrqs; i++) { 142 j = pa[i]; 143 iptr += ict; 144 outdat[j] = iptr; 145 ict = w1[j]; 146 } 147 } 148 149 /* Form the outgoing messages */ 150 /*plug in the headers*/ 151 for (i=0; i<nrqs; i++) { 152 j = pa[i]; 153 outdat[j][0] = 0; 154 ierr = PetscMemzero(outdat[j]+1,2*w3[j]*sizeof(PetscInt));CHKERRQ(ierr); 155 ptr[j] = outdat[j] + 2*w3[j] + 1; 156 } 157 158 /* Memory for doing local proc's work*/ 159 { 160 ierr = PetscCalloc5(imax,&table, imax,&data, imax,&isz, Mbs*imax,&d_p, (Mbs/PETSC_BITS_PER_BYTE+1)*imax,&t_p);CHKERRQ(ierr); 161 162 for (i=0; i<imax; i++) { 163 table[i] = t_p + (Mbs/PETSC_BITS_PER_BYTE+1)*i; 164 data[i] = d_p + (Mbs)*i; 165 } 166 } 167 168 /* Parse the IS and update local tables and the outgoing buf with the data*/ 169 { 170 PetscInt n_i,*data_i,isz_i,*outdat_j,ctr_j; 171 PetscBT table_i; 172 173 for (i=0; i<imax; i++) { 174 ierr = PetscMemzero(ctr,size*sizeof(PetscInt));CHKERRQ(ierr); 175 n_i = n[i]; 176 table_i = table[i]; 177 idx_i = idx[i]; 178 data_i = data[i]; 179 isz_i = isz[i]; 180 for (j=0; j<n_i; j++) { /* parse the indices of each IS */ 181 row = idx_i[j]; 182 ierr = PetscLayoutFindOwner(C->rmap,row*C->rmap->bs,&proc);CHKERRQ(ierr); 183 if (proc != rank) { /* copy to the outgoing buffer */ 184 ctr[proc]++; 185 *ptr[proc] = row; 186 ptr[proc]++; 187 } else { /* Update the local table */ 188 if (!PetscBTLookupSet(table_i,row)) data_i[isz_i++] = row; 189 } 190 } 191 /* Update the headers for the current IS */ 192 for (j=0; j<size; j++) { /* Can Optimise this loop by using pa[] */ 193 if ((ctr_j = ctr[j])) { 194 outdat_j = outdat[j]; 195 k = ++outdat_j[0]; 196 outdat_j[2*k] = ctr_j; 197 outdat_j[2*k-1] = i; 198 } 199 } 200 isz[i] = isz_i; 201 } 202 } 203 204 /* Now post the sends */ 205 ierr = PetscMalloc1(nrqs+1,&s_waits1);CHKERRQ(ierr); 206 for (i=0; i<nrqs; ++i) { 207 j = pa[i]; 208 ierr = MPI_Isend(outdat[j],w1[j],MPIU_INT,j,tag1,comm,s_waits1+i);CHKERRQ(ierr); 209 } 210 211 /* No longer need the original indices*/ 212 for (i=0; i<imax; ++i) { 213 ierr = ISRestoreIndices(is[i],idx+i);CHKERRQ(ierr); 214 } 215 ierr = PetscFree2(idx,n);CHKERRQ(ierr); 216 217 for (i=0; i<imax; ++i) { 218 ierr = ISDestroy(&is[i]);CHKERRQ(ierr); 219 } 220 221 /* Do Local work*/ 222 ierr = MatIncreaseOverlap_MPIBAIJ_Local(C,imax,table,isz,data);CHKERRQ(ierr); 223 224 /* Receive messages*/ 225 ierr = PetscMalloc1(nrqr+1,&recv_status);CHKERRQ(ierr); 226 if (nrqr) {ierr = MPI_Waitall(nrqr,r_waits1,recv_status);CHKERRQ(ierr);} 227 228 ierr = PetscMalloc1(nrqs+1,&s_status);CHKERRQ(ierr); 229 if (nrqs) {ierr = MPI_Waitall(nrqs,s_waits1,s_status);CHKERRQ(ierr);} 230 231 /* Phase 1 sends are complete - deallocate buffers */ 232 ierr = PetscFree4(outdat,ptr,tmp,ctr);CHKERRQ(ierr); 233 ierr = PetscFree4(w1,w2,w3,w4);CHKERRQ(ierr); 234 235 ierr = PetscMalloc1(nrqr+1,&xdata);CHKERRQ(ierr); 236 ierr = PetscMalloc1(nrqr+1,&isz1);CHKERRQ(ierr); 237 ierr = MatIncreaseOverlap_MPIBAIJ_Receive(C,nrqr,rbuf,xdata,isz1);CHKERRQ(ierr); 238 ierr = PetscFree(rbuf[0]);CHKERRQ(ierr); 239 ierr = PetscFree(rbuf);CHKERRQ(ierr); 240 241 /* Send the data back*/ 242 /* Do a global reduction to know the buffer space req for incoming messages*/ 243 { 244 PetscMPIInt *rw1; 245 246 ierr = PetscCalloc1(size,&rw1);CHKERRQ(ierr); 247 248 for (i=0; i<nrqr; ++i) { 249 proc = recv_status[i].MPI_SOURCE; 250 if (proc != onodes1[i]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"MPI_SOURCE mismatch"); 251 rw1[proc] = isz1[i]; 252 } 253 254 ierr = PetscFree(onodes1);CHKERRQ(ierr); 255 ierr = PetscFree(olengths1);CHKERRQ(ierr); 256 257 /* Determine the number of messages to expect, their lengths, from from-ids */ 258 ierr = PetscGatherMessageLengths(comm,nrqr,nrqs,rw1,&onodes2,&olengths2);CHKERRQ(ierr); 259 ierr = PetscFree(rw1);CHKERRQ(ierr); 260 } 261 /* Now post the Irecvs corresponding to these messages */ 262 ierr = PetscPostIrecvInt(comm,tag2,nrqs,onodes2,olengths2,&rbuf2,&r_waits2);CHKERRQ(ierr); 263 264 /* Now post the sends */ 265 ierr = PetscMalloc1(nrqr+1,&s_waits2);CHKERRQ(ierr); 266 for (i=0; i<nrqr; ++i) { 267 j = recv_status[i].MPI_SOURCE; 268 ierr = MPI_Isend(xdata[i],isz1[i],MPIU_INT,j,tag2,comm,s_waits2+i);CHKERRQ(ierr); 269 } 270 271 /* receive work done on other processors*/ 272 { 273 PetscMPIInt idex; 274 PetscInt is_no,ct1,max,*rbuf2_i,isz_i,*data_i,jmax; 275 PetscBT table_i; 276 MPI_Status *status2; 277 278 ierr = PetscMalloc1(PetscMax(nrqr,nrqs)+1,&status2);CHKERRQ(ierr); 279 for (i=0; i<nrqs; ++i) { 280 ierr = MPI_Waitany(nrqs,r_waits2,&idex,status2+i);CHKERRQ(ierr); 281 /* Process the message*/ 282 rbuf2_i = rbuf2[idex]; 283 ct1 = 2*rbuf2_i[0]+1; 284 jmax = rbuf2[idex][0]; 285 for (j=1; j<=jmax; j++) { 286 max = rbuf2_i[2*j]; 287 is_no = rbuf2_i[2*j-1]; 288 isz_i = isz[is_no]; 289 data_i = data[is_no]; 290 table_i = table[is_no]; 291 for (k=0; k<max; k++,ct1++) { 292 row = rbuf2_i[ct1]; 293 if (!PetscBTLookupSet(table_i,row)) data_i[isz_i++] = row; 294 } 295 isz[is_no] = isz_i; 296 } 297 } 298 if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits2,status2);CHKERRQ(ierr);} 299 ierr = PetscFree(status2);CHKERRQ(ierr); 300 } 301 302 for (i=0; i<imax; ++i) { 303 ierr = ISCreateGeneral(PETSC_COMM_SELF,isz[i],data[i],PETSC_COPY_VALUES,is+i);CHKERRQ(ierr); 304 } 305 306 307 ierr = PetscFree(onodes2);CHKERRQ(ierr); 308 ierr = PetscFree(olengths2);CHKERRQ(ierr); 309 310 ierr = PetscFree(pa);CHKERRQ(ierr); 311 ierr = PetscFree(rbuf2[0]);CHKERRQ(ierr); 312 ierr = PetscFree(rbuf2);CHKERRQ(ierr); 313 ierr = PetscFree(s_waits1);CHKERRQ(ierr); 314 ierr = PetscFree(r_waits1);CHKERRQ(ierr); 315 ierr = PetscFree(s_waits2);CHKERRQ(ierr); 316 ierr = PetscFree(r_waits2);CHKERRQ(ierr); 317 ierr = PetscFree5(table,data,isz,d_p,t_p);CHKERRQ(ierr); 318 ierr = PetscFree(s_status);CHKERRQ(ierr); 319 ierr = PetscFree(recv_status);CHKERRQ(ierr); 320 ierr = PetscFree(xdata[0]);CHKERRQ(ierr); 321 ierr = PetscFree(xdata);CHKERRQ(ierr); 322 ierr = PetscFree(isz1);CHKERRQ(ierr); 323 PetscFunctionReturn(0); 324 } 325 326 /* 327 MatIncreaseOverlap_MPIBAIJ_Local - Called by MatincreaseOverlap, to do 328 the work on the local processor. 329 330 Inputs: 331 C - MAT_MPIBAIJ; 332 imax - total no of index sets processed at a time; 333 table - an array of char - size = Mbs bits. 334 335 Output: 336 isz - array containing the count of the solution elements corresponding 337 to each index set; 338 data - pointer to the solutions 339 */ 340 static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Local(Mat C,PetscInt imax,PetscBT *table,PetscInt *isz,PetscInt **data) 341 { 342 Mat_MPIBAIJ *c = (Mat_MPIBAIJ*)C->data; 343 Mat A = c->A,B = c->B; 344 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)B->data; 345 PetscInt start,end,val,max,rstart,cstart,*ai,*aj; 346 PetscInt *bi,*bj,*garray,i,j,k,row,*data_i,isz_i; 347 PetscBT table_i; 348 349 PetscFunctionBegin; 350 rstart = c->rstartbs; 351 cstart = c->cstartbs; 352 ai = a->i; 353 aj = a->j; 354 bi = b->i; 355 bj = b->j; 356 garray = c->garray; 357 358 359 for (i=0; i<imax; i++) { 360 data_i = data[i]; 361 table_i = table[i]; 362 isz_i = isz[i]; 363 for (j=0,max=isz[i]; j<max; j++) { 364 row = data_i[j] - rstart; 365 start = ai[row]; 366 end = ai[row+1]; 367 for (k=start; k<end; k++) { /* Amat */ 368 val = aj[k] + cstart; 369 if (!PetscBTLookupSet(table_i,val)) data_i[isz_i++] = val; 370 } 371 start = bi[row]; 372 end = bi[row+1]; 373 for (k=start; k<end; k++) { /* Bmat */ 374 val = garray[bj[k]]; 375 if (!PetscBTLookupSet(table_i,val)) data_i[isz_i++] = val; 376 } 377 } 378 isz[i] = isz_i; 379 } 380 PetscFunctionReturn(0); 381 } 382 /* 383 MatIncreaseOverlap_MPIBAIJ_Receive - Process the recieved messages, 384 and return the output 385 386 Input: 387 C - the matrix 388 nrqr - no of messages being processed. 389 rbuf - an array of pointers to the recieved requests 390 391 Output: 392 xdata - array of messages to be sent back 393 isz1 - size of each message 394 395 For better efficiency perhaps we should malloc separately each xdata[i], 396 then if a remalloc is required we need only copy the data for that one row 397 rather than all previous rows as it is now where a single large chunck of 398 memory is used. 399 400 */ 401 static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Receive(Mat C,PetscInt nrqr,PetscInt **rbuf,PetscInt **xdata,PetscInt * isz1) 402 { 403 Mat_MPIBAIJ *c = (Mat_MPIBAIJ*)C->data; 404 Mat A = c->A,B = c->B; 405 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)B->data; 406 PetscErrorCode ierr; 407 PetscInt rstart,cstart,*ai,*aj,*bi,*bj,*garray,i,j,k; 408 PetscInt row,total_sz,ct,ct1,ct2,ct3,mem_estimate,oct2,l,start,end; 409 PetscInt val,max1,max2,Mbs,no_malloc =0,*tmp,new_estimate,ctr; 410 PetscInt *rbuf_i,kmax,rbuf_0; 411 PetscBT xtable; 412 413 PetscFunctionBegin; 414 Mbs = c->Mbs; 415 rstart = c->rstartbs; 416 cstart = c->cstartbs; 417 ai = a->i; 418 aj = a->j; 419 bi = b->i; 420 bj = b->j; 421 garray = c->garray; 422 423 424 for (i=0,ct=0,total_sz=0; i<nrqr; ++i) { 425 rbuf_i = rbuf[i]; 426 rbuf_0 = rbuf_i[0]; 427 ct += rbuf_0; 428 for (j=1; j<=rbuf_0; j++) total_sz += rbuf_i[2*j]; 429 } 430 431 if (c->Mbs) max1 = ct*(a->nz +b->nz)/c->Mbs; 432 else max1 = 1; 433 mem_estimate = 3*((total_sz > max1 ? total_sz : max1)+1); 434 ierr = PetscMalloc1(mem_estimate,&xdata[0]);CHKERRQ(ierr); 435 ++no_malloc; 436 ierr = PetscBTCreate(Mbs,&xtable);CHKERRQ(ierr); 437 ierr = PetscMemzero(isz1,nrqr*sizeof(PetscInt));CHKERRQ(ierr); 438 439 ct3 = 0; 440 for (i=0; i<nrqr; i++) { /* for easch mesg from proc i */ 441 rbuf_i = rbuf[i]; 442 rbuf_0 = rbuf_i[0]; 443 ct1 = 2*rbuf_0+1; 444 ct2 = ct1; 445 ct3 += ct1; 446 for (j=1; j<=rbuf_0; j++) { /* for each IS from proc i*/ 447 ierr = PetscBTMemzero(Mbs,xtable);CHKERRQ(ierr); 448 oct2 = ct2; 449 kmax = rbuf_i[2*j]; 450 for (k=0; k<kmax; k++,ct1++) { 451 row = rbuf_i[ct1]; 452 if (!PetscBTLookupSet(xtable,row)) { 453 if (!(ct3 < mem_estimate)) { 454 new_estimate = (PetscInt)(1.5*mem_estimate)+1; 455 ierr = PetscMalloc1(new_estimate,&tmp);CHKERRQ(ierr); 456 ierr = PetscMemcpy(tmp,xdata[0],mem_estimate*sizeof(PetscInt));CHKERRQ(ierr); 457 ierr = PetscFree(xdata[0]);CHKERRQ(ierr); 458 xdata[0] = tmp; 459 mem_estimate = new_estimate; ++no_malloc; 460 for (ctr=1; ctr<=i; ctr++) xdata[ctr] = xdata[ctr-1] + isz1[ctr-1]; 461 } 462 xdata[i][ct2++] = row; 463 ct3++; 464 } 465 } 466 for (k=oct2,max2=ct2; k<max2; k++) { 467 row = xdata[i][k] - rstart; 468 start = ai[row]; 469 end = ai[row+1]; 470 for (l=start; l<end; l++) { 471 val = aj[l] + cstart; 472 if (!PetscBTLookupSet(xtable,val)) { 473 if (!(ct3 < mem_estimate)) { 474 new_estimate = (PetscInt)(1.5*mem_estimate)+1; 475 ierr = PetscMalloc1(new_estimate,&tmp);CHKERRQ(ierr); 476 ierr = PetscMemcpy(tmp,xdata[0],mem_estimate*sizeof(PetscInt));CHKERRQ(ierr); 477 ierr = PetscFree(xdata[0]);CHKERRQ(ierr); 478 xdata[0] = tmp; 479 mem_estimate = new_estimate; ++no_malloc; 480 for (ctr=1; ctr<=i; ctr++) xdata[ctr] = xdata[ctr-1] + isz1[ctr-1]; 481 } 482 xdata[i][ct2++] = val; 483 ct3++; 484 } 485 } 486 start = bi[row]; 487 end = bi[row+1]; 488 for (l=start; l<end; l++) { 489 val = garray[bj[l]]; 490 if (!PetscBTLookupSet(xtable,val)) { 491 if (!(ct3 < mem_estimate)) { 492 new_estimate = (PetscInt)(1.5*mem_estimate)+1; 493 ierr = PetscMalloc1(new_estimate,&tmp);CHKERRQ(ierr); 494 ierr = PetscMemcpy(tmp,xdata[0],mem_estimate*sizeof(PetscInt));CHKERRQ(ierr); 495 ierr = PetscFree(xdata[0]);CHKERRQ(ierr); 496 xdata[0] = tmp; 497 mem_estimate = new_estimate; ++no_malloc; 498 for (ctr =1; ctr <=i; ctr++) xdata[ctr] = xdata[ctr-1] + isz1[ctr-1]; 499 } 500 xdata[i][ct2++] = val; 501 ct3++; 502 } 503 } 504 } 505 /* Update the header*/ 506 xdata[i][2*j] = ct2 - oct2; /* Undo the vector isz1 and use only a var*/ 507 xdata[i][2*j-1] = rbuf_i[2*j-1]; 508 } 509 xdata[i][0] = rbuf_0; 510 xdata[i+1] = xdata[i] + ct2; 511 isz1[i] = ct2; /* size of each message */ 512 } 513 ierr = PetscBTDestroy(&xtable);CHKERRQ(ierr); 514 ierr = PetscInfo3(C,"Allocated %D bytes, required %D, no of mallocs = %D\n",mem_estimate,ct3,no_malloc);CHKERRQ(ierr); 515 PetscFunctionReturn(0); 516 } 517 518 PetscErrorCode MatGetSubMatrices_MPIBAIJ(Mat C,PetscInt ismax,const IS isrow[],const IS iscol[],MatReuse scall,Mat *submat[]) 519 { 520 IS *isrow_block,*iscol_block; 521 Mat_MPIBAIJ *c = (Mat_MPIBAIJ*)C->data; 522 PetscErrorCode ierr; 523 PetscInt nmax,nstages_local,nstages,i,pos,max_no,ncol,nrow,N=C->cmap->N,bs=C->rmap->bs; 524 PetscBool colflag,*allcolumns,*allrows; 525 526 PetscFunctionBegin; 527 //printf("MatGetSubMatrices_MPIBAIJ scall %d, bs %d\n",scall,bs); 528 /* The compression and expansion should be avoided. Doesn't point 529 out errors, might change the indices, hence buggey */ 530 ierr = PetscMalloc2(ismax+1,&isrow_block,ismax+1,&iscol_block);CHKERRQ(ierr); 531 ierr = ISCompressIndicesGeneral(N,C->rmap->n,bs,ismax,isrow,isrow_block);CHKERRQ(ierr); 532 ierr = ISCompressIndicesGeneral(N,C->cmap->n,bs,ismax,iscol,iscol_block);CHKERRQ(ierr); 533 534 /* Check for special case: each processor gets entire matrix columns -- rm! */ 535 ierr = PetscMalloc2(ismax+1,&allcolumns,ismax+1,&allrows);CHKERRQ(ierr); 536 for (i=0; i<ismax; i++) { 537 ierr = ISIdentity(iscol[i],&colflag);CHKERRQ(ierr); 538 ierr = ISGetLocalSize(iscol[i],&ncol);CHKERRQ(ierr); 539 if (colflag && ncol == C->cmap->N) allcolumns[i] = PETSC_TRUE; 540 else allcolumns[i] = PETSC_FALSE; 541 542 ierr = ISIdentity(isrow[i],&colflag);CHKERRQ(ierr); 543 ierr = ISGetLocalSize(isrow[i],&nrow);CHKERRQ(ierr); 544 if (colflag && nrow == C->rmap->N) allrows[i] = PETSC_TRUE; 545 else allrows[i] = PETSC_FALSE; 546 } 547 548 /* Allocate memory to hold all the submatrices */ 549 if (scall == MAT_INITIAL_MATRIX) { 550 ierr = PetscMalloc1(ismax+1,submat);CHKERRQ(ierr); 551 } 552 /* Determine the number of stages through which submatrices are done */ 553 nmax = 20*1000000 / (c->Nbs * sizeof(PetscInt)); 554 if (!nmax) nmax = 1; 555 nstages_local = ismax/nmax + ((ismax % nmax) ? 1 : 0); 556 557 /* Make sure every processor loops through the nstages */ 558 ierr = MPIU_Allreduce(&nstages_local,&nstages,1,MPIU_INT,MPI_MAX,PetscObjectComm((PetscObject)C));CHKERRQ(ierr); 559 for (i=0,pos=0; i<nstages; i++) { 560 if (pos+nmax <= ismax) max_no = nmax; 561 else if (pos == ismax) max_no = 0; 562 else max_no = ismax-pos; 563 564 PetscBool isbaij; 565 ierr = PetscObjectTypeCompare((PetscObject)C,MATMPIBAIJ,&isbaij);CHKERRQ(ierr); 566 if (isbaij) { 567 ierr = MatGetSubMatrices_MPIBAIJ_local_new(C,max_no,isrow_block+pos,iscol_block+pos,scall,allrows+pos,allcolumns+pos,*submat+pos);CHKERRQ(ierr); 568 } else { 569 ierr = MatGetSubMatrices_MPIBAIJ_local(C,max_no,isrow_block+pos,iscol_block+pos,scall,allrows+pos,allcolumns+pos,*submat+pos);CHKERRQ(ierr); 570 } 571 pos += max_no; 572 } 573 574 for (i=0; i<ismax; i++) { 575 ierr = ISDestroy(&isrow_block[i]);CHKERRQ(ierr); 576 ierr = ISDestroy(&iscol_block[i]);CHKERRQ(ierr); 577 } 578 ierr = PetscFree2(isrow_block,iscol_block);CHKERRQ(ierr); 579 ierr = PetscFree2(allcolumns,allrows);CHKERRQ(ierr); 580 PetscFunctionReturn(0); 581 } 582 583 #if defined(PETSC_USE_CTABLE) 584 PetscErrorCode PetscGetProc(const PetscInt row, const PetscMPIInt size, const PetscInt proc_gnode[], PetscMPIInt *rank) 585 { 586 PetscInt nGlobalNd = proc_gnode[size]; 587 PetscMPIInt fproc; 588 PetscErrorCode ierr; 589 590 PetscFunctionBegin; 591 ierr = PetscMPIIntCast((PetscInt)(((float)row * (float)size / (float)nGlobalNd + 0.5)),&fproc);CHKERRQ(ierr); 592 if (fproc > size) fproc = size; 593 while (row < proc_gnode[fproc] || row >= proc_gnode[fproc+1]) { 594 if (row < proc_gnode[fproc]) fproc--; 595 else fproc++; 596 } 597 *rank = fproc; 598 PetscFunctionReturn(0); 599 } 600 #endif 601 602 /* -------------------------------------------------------------------------*/ 603 /* This code is used for BAIJ and SBAIJ matrices (unfortunate dependency) */ 604 605 PetscErrorCode MatDestroy_MPIBAIJ_MatGetSubmatrices(Mat C) 606 { 607 PetscErrorCode ierr; 608 Mat_SeqBAIJ *c = (Mat_SeqBAIJ*)C->data; 609 Mat_SubMat *submatj = c->submatis1; 610 PetscInt i; 611 612 PetscFunctionBegin; 613 if (!submatj->id) { /* delete data that are linked only to submats[id=0] */ 614 ierr = PetscFree4(submatj->sbuf1,submatj->ptr,submatj->tmp,submatj->ctr);CHKERRQ(ierr); 615 616 for (i=0; i<submatj->nrqr; ++i) { 617 ierr = PetscFree(submatj->sbuf2[i]);CHKERRQ(ierr); 618 } 619 ierr = PetscFree3(submatj->sbuf2,submatj->req_size,submatj->req_source1);CHKERRQ(ierr); 620 621 if (submatj->rbuf1) { 622 ierr = PetscFree(submatj->rbuf1[0]);CHKERRQ(ierr); 623 ierr = PetscFree(submatj->rbuf1);CHKERRQ(ierr); 624 } 625 626 for (i=0; i<submatj->nrqs; ++i) { 627 ierr = PetscFree(submatj->rbuf3[i]);CHKERRQ(ierr); 628 } 629 ierr = PetscFree3(submatj->req_source2,submatj->rbuf2,submatj->rbuf3);CHKERRQ(ierr); 630 ierr = PetscFree(submatj->pa);CHKERRQ(ierr); 631 } 632 633 #if defined(PETSC_USE_CTABLE) 634 ierr = PetscTableDestroy((PetscTable*)&submatj->rmap);CHKERRQ(ierr); 635 if (submatj->cmap_loc) {ierr = PetscFree(submatj->cmap_loc);CHKERRQ(ierr);} 636 ierr = PetscFree(submatj->rmap_loc);CHKERRQ(ierr); 637 #else 638 ierr = PetscFree(submatj->rmap);CHKERRQ(ierr); 639 #endif 640 641 if (!submatj->allcolumns) { 642 #if defined(PETSC_USE_CTABLE) 643 ierr = PetscTableDestroy((PetscTable*)&submatj->cmap);CHKERRQ(ierr); 644 #else 645 ierr = PetscFree(submatj->cmap);CHKERRQ(ierr); 646 #endif 647 } 648 ierr = submatj->destroy(C);CHKERRQ(ierr); 649 ierr = PetscFree(submatj->row2proc);CHKERRQ(ierr); 650 651 ierr = PetscFree(submatj);CHKERRQ(ierr); 652 PetscFunctionReturn(0); 653 } 654 655 PetscErrorCode MatGetSubMatrices_MPIBAIJ_local_new(Mat C,PetscInt ismax,const IS isrow[],const IS iscol[],MatReuse scall,PetscBool *allrows,PetscBool *allcolumns,Mat *submats) 656 { 657 Mat_MPIBAIJ *c = (Mat_MPIBAIJ*)C->data; 658 Mat A = c->A; 659 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)c->B->data,*subc; 660 const PetscInt **icol,**irow; 661 PetscInt *nrow,*ncol,start; 662 PetscErrorCode ierr; 663 PetscMPIInt rank,size,tag0,tag2,tag3,tag4,*w1,*w2,*w3,*w4,nrqr; 664 PetscInt **sbuf1,**sbuf2,i,j,k,l,ct1,ct2,**rbuf1,row,proc=-1; 665 PetscInt nrqs=0,msz,**ptr=NULL,*req_size=NULL,*ctr=NULL,*pa,*tmp=NULL,tcol; 666 PetscInt **rbuf3=NULL,*req_source1=NULL,*req_source2,**sbuf_aj,**rbuf2=NULL,max1,max2; 667 PetscInt **lens,is_no,ncols,*cols,mat_i,*mat_j,tmp2,jmax; 668 #if defined(PETSC_USE_CTABLE) 669 PetscTable *cmap,cmap_i=NULL,*rmap,rmap_i; 670 #else 671 PetscInt **cmap,*cmap_i=NULL,**rmap,*rmap_i; 672 #endif 673 const PetscInt *irow_i; 674 PetscInt ctr_j,*sbuf1_j,*sbuf_aj_i,*rbuf1_i,kmax,*lens_i; 675 MPI_Request *s_waits1,*r_waits1,*s_waits2,*r_waits2,*r_waits3; 676 MPI_Request *r_waits4,*s_waits3,*s_waits4; 677 MPI_Status *r_status1,*r_status2,*s_status1,*s_status3,*s_status2; 678 MPI_Status *r_status3,*r_status4,*s_status4; 679 MPI_Comm comm; 680 PetscScalar **rbuf4,*rbuf4_i,**sbuf_aa,*vals,*mat_a,*imat_a,*sbuf_aa_i; 681 PetscMPIInt *onodes1,*olengths1,end; 682 PetscInt **row2proc,*row2proc_i,ilen_row,*imat_ilen,*imat_j,*imat_i; 683 Mat_SubMat **smats,*smat_i; 684 PetscBool *issorted,colflag,iscsorted=PETSC_TRUE; 685 PetscInt *sbuf1_i,*rbuf2_i,*rbuf3_i,ilen; 686 687 PetscInt bs=C->rmap->bs,bs2=c->bs2,rstart = c->rstartbs; 688 PetscBool ijonly=c->ijonly; /* private flag indicates only matrix data structures are requested */ 689 PetscInt nzA,nzB,*a_i=a->i,*b_i=b->i,*a_j = a->j,*b_j = b->j,ctmp,imark,*cworkA,*cworkB; 690 PetscScalar *vworkA,*vworkB,*a_a = a->a,*b_a = b->a; 691 PetscInt cstart = c->cstartbs,*bmap = c->garray; 692 693 PetscFunctionBegin; 694 ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr); 695 size = c->size; 696 rank = c->rank; 697 698 ierr = PetscMalloc5(ismax,&smats,ismax,&row2proc,ismax,&cmap,ismax,&rmap,ismax,&allcolumns);CHKERRQ(ierr); 699 ierr = PetscMalloc5(ismax,&irow,ismax,&icol,ismax,&nrow,ismax,&ncol,ismax,&issorted);CHKERRQ(ierr); 700 701 for (i=0; i<ismax; i++) { 702 ierr = ISSorted(iscol[i],&issorted[i]);CHKERRQ(ierr); 703 if (!issorted[i]) iscsorted = issorted[i]; 704 705 ierr = ISSorted(isrow[i],&issorted[i]);CHKERRQ(ierr); 706 707 ierr = ISGetIndices(isrow[i],&irow[i]);CHKERRQ(ierr); 708 ierr = ISGetLocalSize(isrow[i],&nrow[i]);CHKERRQ(ierr); 709 710 /* Check for special case: allcolumn */ 711 ierr = ISIdentity(iscol[i],&colflag);CHKERRQ(ierr); 712 ierr = ISGetLocalSize(iscol[i],&ncol[i]);CHKERRQ(ierr); 713 714 if (colflag && ncol[i] == c->Nbs) { 715 allcolumns[i] = PETSC_TRUE; 716 icol[i] = NULL; 717 printf("[%d] allcolumns[%d] true\n",rank,i); 718 } else { 719 allcolumns[i] = PETSC_FALSE; 720 ierr = ISGetIndices(iscol[i],&icol[i]);CHKERRQ(ierr); 721 } 722 } 723 //printf("[%d] nrow %d, ncol %d\n",rank,nrow[0],ncol[0]); 724 725 if (scall == MAT_REUSE_MATRIX) { 726 /* Assumes new rows are same length as the old rows */ 727 for (i=0; i<ismax; i++) { 728 subc = (Mat_SeqBAIJ*)(submats[i]->data); 729 if (subc->mbs != nrow[i] || subc->nbs != ncol[i]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong size"); 730 731 /* Initial matrix as if empty */ 732 ierr = PetscMemzero(subc->ilen,subc->mbs*sizeof(PetscInt));CHKERRQ(ierr); 733 734 /* Initial matrix as if empty */ 735 submats[i]->factortype = C->factortype; 736 737 smat_i = subc->submatis1; 738 smats[i] = smat_i; 739 740 nrqs = smat_i->nrqs; 741 nrqr = smat_i->nrqr; 742 rbuf1 = smat_i->rbuf1; 743 rbuf2 = smat_i->rbuf2; 744 rbuf3 = smat_i->rbuf3; 745 req_source2 = smat_i->req_source2; 746 747 sbuf1 = smat_i->sbuf1; 748 sbuf2 = smat_i->sbuf2; 749 ptr = smat_i->ptr; 750 tmp = smat_i->tmp; 751 ctr = smat_i->ctr; 752 753 pa = smat_i->pa; 754 req_size = smat_i->req_size; 755 req_source1 = smat_i->req_source1; 756 757 allcolumns[i] = smat_i->allcolumns; 758 row2proc[i] = smat_i->row2proc; 759 rmap[i] = smat_i->rmap; 760 cmap[i] = smat_i->cmap; 761 } 762 } else { /* scall == MAT_INITIAL_MATRIX */ 763 /* Get some new tags to keep the communication clean */ 764 ierr = PetscObjectGetNewTag((PetscObject)C,&tag2);CHKERRQ(ierr); 765 ierr = PetscObjectGetNewTag((PetscObject)C,&tag3);CHKERRQ(ierr); 766 767 /* evaluate communication - mesg to who, length of mesg, and buffer space 768 required. Based on this, buffers are allocated, and data copied into them*/ 769 ierr = PetscCalloc4(size,&w1,size,&w2,size,&w3,size,&w4);CHKERRQ(ierr); /* mesg size, initialize work vectors */ 770 771 for (i=0; i<ismax; i++) { 772 jmax = nrow[i]; 773 irow_i = irow[i]; 774 775 ierr = PetscMalloc1(jmax,&row2proc_i);CHKERRQ(ierr); 776 row2proc[i] = row2proc_i; 777 778 if (issorted[i]) proc = 0; 779 for (j=0; j<jmax; j++) { 780 if (!issorted[i]) proc = 0; 781 row = irow_i[j]; 782 while (row >= c->rangebs[proc+1]) proc++; 783 w4[proc]++; 784 row2proc_i[j] = proc; /* map row index to proc */ 785 } 786 for (j=0; j<size; j++) { 787 if (w4[j]) { w1[j] += w4[j]; w3[j]++; w4[j] = 0;} 788 } 789 } 790 791 nrqs = 0; /* no of outgoing messages */ 792 msz = 0; /* total mesg length (for all procs) */ 793 w1[rank] = 0; /* no mesg sent to self */ 794 w3[rank] = 0; 795 for (i=0; i<size; i++) { 796 if (w1[i]) { w2[i] = 1; nrqs++;} /* there exists a message to proc i */ 797 } 798 ierr = PetscMalloc1(nrqs+1,&pa);CHKERRQ(ierr); /*(proc -array)*/ 799 for (i=0,j=0; i<size; i++) { 800 if (w1[i]) { pa[j] = i; j++; } 801 } 802 803 /* Each message would have a header = 1 + 2*(no of IS) + data */ 804 for (i=0; i<nrqs; i++) { 805 j = pa[i]; 806 w1[j] += w2[j] + 2* w3[j]; 807 msz += w1[j]; 808 } 809 ierr = PetscInfo2(0,"Number of outgoing messages %D Total message length %D\n",nrqs,msz);CHKERRQ(ierr); 810 811 /* Determine the number of messages to expect, their lengths, from from-ids */ 812 ierr = PetscGatherNumberOfMessages(comm,w2,w1,&nrqr);CHKERRQ(ierr); 813 ierr = PetscGatherMessageLengths(comm,nrqs,nrqr,w1,&onodes1,&olengths1);CHKERRQ(ierr); 814 815 /* Now post the Irecvs corresponding to these messages */ 816 tag0 = ((PetscObject)C)->tag; 817 ierr = PetscPostIrecvInt(comm,tag0,nrqr,onodes1,olengths1,&rbuf1,&r_waits1);CHKERRQ(ierr); 818 /* printf("[%d] nrqs %d, nrqr %d\n",rank,nrqs,nrqr); */ 819 820 ierr = PetscFree(onodes1);CHKERRQ(ierr); 821 ierr = PetscFree(olengths1);CHKERRQ(ierr); 822 823 /* Allocate Memory for outgoing messages */ 824 ierr = PetscMalloc4(size,&sbuf1,size,&ptr,2*msz,&tmp,size,&ctr);CHKERRQ(ierr); 825 ierr = PetscMemzero(sbuf1,size*sizeof(PetscInt*));CHKERRQ(ierr); 826 ierr = PetscMemzero(ptr,size*sizeof(PetscInt*));CHKERRQ(ierr); 827 828 { 829 PetscInt *iptr = tmp; 830 k = 0; 831 for (i=0; i<nrqs; i++) { 832 j = pa[i]; 833 iptr += k; 834 sbuf1[j] = iptr; 835 k = w1[j]; 836 } 837 } 838 839 /* Form the outgoing messages. Initialize the header space */ 840 for (i=0; i<nrqs; i++) { 841 j = pa[i]; 842 sbuf1[j][0] = 0; 843 ierr = PetscMemzero(sbuf1[j]+1,2*w3[j]*sizeof(PetscInt));CHKERRQ(ierr); 844 ptr[j] = sbuf1[j] + 2*w3[j] + 1; 845 } 846 847 /* Parse the isrow and copy data into outbuf */ 848 for (i=0; i<ismax; i++) { 849 row2proc_i = row2proc[i]; 850 ierr = PetscMemzero(ctr,size*sizeof(PetscInt));CHKERRQ(ierr); 851 irow_i = irow[i]; 852 jmax = nrow[i]; 853 for (j=0; j<jmax; j++) { /* parse the indices of each IS */ 854 proc = row2proc_i[j]; 855 if (proc != rank) { /* copy to the outgoing buf*/ 856 ctr[proc]++; 857 *ptr[proc] = irow_i[j]; 858 ptr[proc]++; 859 } 860 } 861 /* Update the headers for the current IS */ 862 for (j=0; j<size; j++) { /* Can Optimise this loop too */ 863 if ((ctr_j = ctr[j])) { 864 sbuf1_j = sbuf1[j]; 865 k = ++sbuf1_j[0]; 866 sbuf1_j[2*k] = ctr_j; 867 sbuf1_j[2*k-1] = i; 868 } 869 } 870 } 871 872 /* Now post the sends */ 873 ierr = PetscMalloc1(nrqs+1,&s_waits1);CHKERRQ(ierr); 874 for (i=0; i<nrqs; ++i) { 875 j = pa[i]; 876 ierr = MPI_Isend(sbuf1[j],w1[j],MPIU_INT,j,tag0,comm,s_waits1+i);CHKERRQ(ierr); 877 } 878 879 /* Post Receives to capture the buffer size */ 880 ierr = PetscMalloc1(nrqs+1,&r_waits2);CHKERRQ(ierr); 881 ierr = PetscMalloc3(nrqs+1,&req_source2,nrqs+1,&rbuf2,nrqs+1,&rbuf3);CHKERRQ(ierr); 882 rbuf2[0] = tmp + msz; 883 for (i=1; i<nrqs; ++i) { 884 rbuf2[i] = rbuf2[i-1]+w1[pa[i-1]]; 885 } 886 for (i=0; i<nrqs; ++i) { 887 j = pa[i]; 888 ierr = MPI_Irecv(rbuf2[i],w1[j],MPIU_INT,j,tag2,comm,r_waits2+i);CHKERRQ(ierr); 889 } 890 891 /* Send to other procs the buf size they should allocate */ 892 /* Receive messages*/ 893 ierr = PetscMalloc1(nrqr+1,&s_waits2);CHKERRQ(ierr); 894 ierr = PetscMalloc1(nrqr+1,&r_status1);CHKERRQ(ierr); 895 ierr = PetscMalloc3(nrqr,&sbuf2,nrqr,&req_size,nrqr,&req_source1);CHKERRQ(ierr); 896 { 897 PetscInt *sAi = a->i,*sBi = b->i,id; 898 PetscInt *sbuf2_i; 899 900 ierr = MPI_Waitall(nrqr,r_waits1,r_status1);CHKERRQ(ierr); 901 for (i=0; i<nrqr; ++i) { 902 req_size[i] = 0; 903 rbuf1_i = rbuf1[i]; 904 start = 2*rbuf1_i[0] + 1; 905 ierr = MPI_Get_count(r_status1+i,MPIU_INT,&end);CHKERRQ(ierr); 906 ierr = PetscMalloc1(end+1,&sbuf2[i]);CHKERRQ(ierr); 907 sbuf2_i = sbuf2[i]; 908 for (j=start; j<end; j++) { 909 id = rbuf1_i[j] - rstart; 910 ncols = sAi[id+1] - sAi[id] + sBi[id+1] - sBi[id]; 911 sbuf2_i[j] = ncols; 912 req_size[i] += ncols; 913 } 914 req_source1[i] = r_status1[i].MPI_SOURCE; 915 /* form the header */ 916 sbuf2_i[0] = req_size[i]; 917 for (j=1; j<start; j++) sbuf2_i[j] = rbuf1_i[j]; 918 919 ierr = MPI_Isend(sbuf2_i,end,MPIU_INT,req_source1[i],tag2,comm,s_waits2+i);CHKERRQ(ierr); 920 } 921 } 922 ierr = PetscFree(r_status1);CHKERRQ(ierr); 923 ierr = PetscFree(r_waits1);CHKERRQ(ierr); 924 ierr = PetscFree4(w1,w2,w3,w4);CHKERRQ(ierr); 925 926 /* Receive messages*/ 927 ierr = PetscMalloc1(nrqs+1,&r_waits3);CHKERRQ(ierr); 928 ierr = PetscMalloc1(nrqs+1,&r_status2);CHKERRQ(ierr); 929 930 ierr = MPI_Waitall(nrqs,r_waits2,r_status2);CHKERRQ(ierr); 931 for (i=0; i<nrqs; ++i) { 932 ierr = PetscMalloc1(rbuf2[i][0]+1,&rbuf3[i]);CHKERRQ(ierr); 933 req_source2[i] = r_status2[i].MPI_SOURCE; 934 ierr = MPI_Irecv(rbuf3[i],rbuf2[i][0],MPIU_INT,req_source2[i],tag3,comm,r_waits3+i);CHKERRQ(ierr); 935 } 936 ierr = PetscFree(r_status2);CHKERRQ(ierr); 937 ierr = PetscFree(r_waits2);CHKERRQ(ierr); 938 939 /* Wait on sends1 and sends2 */ 940 ierr = PetscMalloc1(nrqs+1,&s_status1);CHKERRQ(ierr); 941 ierr = PetscMalloc1(nrqr+1,&s_status2);CHKERRQ(ierr); 942 943 if (nrqs) {ierr = MPI_Waitall(nrqs,s_waits1,s_status1);CHKERRQ(ierr);} 944 if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits2,s_status2);CHKERRQ(ierr);} 945 ierr = PetscFree(s_status1);CHKERRQ(ierr); 946 ierr = PetscFree(s_status2);CHKERRQ(ierr); 947 ierr = PetscFree(s_waits1);CHKERRQ(ierr); 948 ierr = PetscFree(s_waits2);CHKERRQ(ierr); 949 950 /* Now allocate sending buffers for a->j, and send them off */ 951 ierr = PetscMalloc1(nrqr+1,&sbuf_aj);CHKERRQ(ierr); 952 for (i=0,j=0; i<nrqr; i++) j += req_size[i]; 953 ierr = PetscMalloc1(j+1,&sbuf_aj[0]);CHKERRQ(ierr); 954 for (i=1; i<nrqr; i++) sbuf_aj[i] = sbuf_aj[i-1] + req_size[i-1]; 955 956 ierr = PetscMalloc1(nrqr+1,&s_waits3);CHKERRQ(ierr); 957 { 958 959 for (i=0; i<nrqr; i++) { 960 rbuf1_i = rbuf1[i]; 961 sbuf_aj_i = sbuf_aj[i]; 962 ct1 = 2*rbuf1_i[0] + 1; 963 ct2 = 0; 964 for (j=1,max1=rbuf1_i[0]; j<=max1; j++) { 965 kmax = rbuf1[i][2*j]; 966 for (k=0; k<kmax; k++,ct1++) { 967 row = rbuf1_i[ct1] - rstart; 968 nzA = a_i[row+1] - a_i[row]; nzB = b_i[row+1] - b_i[row]; 969 ncols = nzA + nzB; 970 cworkA = a_j + a_i[row]; cworkB = b_j + b_i[row]; 971 972 /* load the column indices for this row into cols */ 973 cols = sbuf_aj_i + ct2; 974 for (l=0; l<nzB; l++) { 975 if ((ctmp = bmap[cworkB[l]]) < cstart) cols[l] = ctmp; 976 else break; 977 } 978 imark = l; 979 for (l=0; l<nzA; l++) {cols[imark+l] = cstart + cworkA[l];} 980 for (l=imark; l<nzB; l++) cols[nzA+l] = bmap[cworkB[l]]; 981 ct2 += ncols; 982 } 983 } 984 ierr = MPI_Isend(sbuf_aj_i,req_size[i],MPIU_INT,req_source1[i],tag3,comm,s_waits3+i);CHKERRQ(ierr); 985 } 986 } 987 ierr = PetscMalloc2(nrqs+1,&r_status3,nrqr+1,&s_status3);CHKERRQ(ierr); 988 989 /* create col map: global col of C -> local col of submatrices */ 990 const PetscInt *icol_i; 991 #if defined(PETSC_USE_CTABLE) 992 for (i=0; i<ismax; i++) { 993 if (!allcolumns[i]) { 994 ierr = PetscTableCreate(ncol[i]+1,c->Nbs+1,&cmap[i]);CHKERRQ(ierr); 995 996 jmax = ncol[i]; 997 icol_i = icol[i]; 998 cmap_i = cmap[i]; 999 for (j=0; j<jmax; j++) { 1000 ierr = PetscTableAdd(cmap[i],icol_i[j]+1,j+1,INSERT_VALUES);CHKERRQ(ierr); 1001 } 1002 } else cmap[i] = NULL; 1003 } 1004 #else 1005 for (i=0; i<ismax; i++) { 1006 if (!allcolumns[i]) { 1007 ierr = PetscCalloc1(c->Nbs,&cmap[i]);CHKERRQ(ierr); 1008 jmax = ncol[i]; 1009 icol_i = icol[i]; 1010 cmap_i = cmap[i]; 1011 for (j=0; j<jmax; j++) cmap_i[icol_i[j]] = j+1; 1012 } else cmap[i] = NULL; 1013 } 1014 #endif 1015 1016 /* Create lens which is required for MatCreate... */ 1017 for (i=0,j=0; i<ismax; i++) j += nrow[i]; 1018 ierr = PetscMalloc1(ismax,&lens);CHKERRQ(ierr); 1019 1020 if (ismax) { 1021 ierr = PetscCalloc1(j,&lens[0]);CHKERRQ(ierr); 1022 } 1023 for (i=1; i<ismax; i++) lens[i] = lens[i-1] + nrow[i-1]; 1024 1025 /* Update lens from local data */ 1026 for (i=0; i<ismax; i++) { 1027 row2proc_i = row2proc[i]; 1028 jmax = nrow[i]; 1029 if (!allcolumns[i]) cmap_i = cmap[i]; 1030 irow_i = irow[i]; 1031 lens_i = lens[i]; 1032 for (j=0; j<jmax; j++) { 1033 row = irow_i[j]; /* global blocked row of C */ 1034 proc = row2proc_i[j]; 1035 if (proc == rank) { 1036 /* Get indices from matA and then from matB */ 1037 #if defined(PETSC_USE_CTABLE) 1038 PetscInt tt; 1039 #endif 1040 row = row - rstart; 1041 nzA = a_i[row+1] - a_i[row]; 1042 nzB = b_i[row+1] - b_i[row]; 1043 cworkA = a_j + a_i[row]; 1044 cworkB = b_j + b_i[row]; 1045 1046 if (!allcolumns[i]) { 1047 #if defined(PETSC_USE_CTABLE) 1048 for (k=0; k<nzA; k++) { 1049 ierr = PetscTableFind(cmap_i,cstart+cworkA[k]+1,&tt);CHKERRQ(ierr); 1050 if (tt) lens_i[j]++; 1051 } 1052 for (k=0; k<nzB; k++) { 1053 ierr = PetscTableFind(cmap_i,bmap[cworkB[k]]+1,&tt);CHKERRQ(ierr); 1054 if (tt) lens_i[j]++; 1055 } 1056 1057 #else 1058 for (k=0; k<nzA; k++) { 1059 if (cmap_i[cstart + cworkA[k]]) lens_i[j]++; 1060 } 1061 for (k=0; k<nzB; k++) { 1062 if (cmap_i[bmap[cworkB[k]]]) lens_i[j]++; 1063 } 1064 #endif 1065 } else { /* allcolumns */ 1066 lens_i[j] = nzA + nzB; 1067 } 1068 } 1069 } 1070 } 1071 1072 /* Create row map: global row of C -> local row of submatrices */ 1073 #if defined(PETSC_USE_CTABLE) 1074 for (i=0; i<ismax; i++) { 1075 ierr = PetscTableCreate(nrow[i]+1,c->Mbs+1,&rmap[i]);CHKERRQ(ierr); 1076 irow_i = irow[i]; 1077 jmax = nrow[i]; 1078 for (j=0; j<jmax; j++) { 1079 ierr = PetscTableAdd(rmap[i],irow_i[j]+1,j+1,INSERT_VALUES);CHKERRQ(ierr); 1080 } 1081 } 1082 #else 1083 for (i=0; i<ismax; i++) { 1084 ierr = PetscCalloc1(c->Mbs,&rmap[i]);CHKERRQ(ierr); 1085 rmap_i = rmap[i]; 1086 irow_i = irow[i]; 1087 jmax = nrow[i]; 1088 for (j=0; j<jmax; j++) { 1089 rmap_i[irow_i[j]] = j; 1090 } 1091 } 1092 #endif 1093 1094 /* Update lens from offproc data */ 1095 { 1096 PetscInt *rbuf2_i,*rbuf3_i,*sbuf1_i; 1097 1098 ierr = MPI_Waitall(nrqs,r_waits3,r_status3);CHKERRQ(ierr); 1099 for (tmp2=0; tmp2<nrqs; tmp2++) { 1100 sbuf1_i = sbuf1[pa[tmp2]]; 1101 jmax = sbuf1_i[0]; 1102 ct1 = 2*jmax+1; 1103 ct2 = 0; 1104 rbuf2_i = rbuf2[tmp2]; 1105 rbuf3_i = rbuf3[tmp2]; 1106 for (j=1; j<=jmax; j++) { 1107 is_no = sbuf1_i[2*j-1]; 1108 max1 = sbuf1_i[2*j]; 1109 lens_i = lens[is_no]; 1110 if (!allcolumns[is_no]) cmap_i = cmap[is_no]; 1111 rmap_i = rmap[is_no]; 1112 for (k=0; k<max1; k++,ct1++) { 1113 #if defined(PETSC_USE_CTABLE) 1114 ierr = PetscTableFind(rmap_i,sbuf1_i[ct1]+1,&row);CHKERRQ(ierr); 1115 row--; 1116 if (row < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"row not found in table"); 1117 #else 1118 row = rmap_i[sbuf1_i[ct1]]; /* the val in the new matrix to be */ 1119 #endif 1120 max2 = rbuf2_i[ct1]; 1121 for (l=0; l<max2; l++,ct2++) { 1122 if (!allcolumns[is_no]) { 1123 #if defined(PETSC_USE_CTABLE) 1124 ierr = PetscTableFind(cmap_i,rbuf3_i[ct2]+1,&tcol);CHKERRQ(ierr); 1125 #else 1126 tcol = cmap_i[rbuf3_i[ct2]]; 1127 #endif 1128 if (tcol) lens_i[row]++; 1129 } else { /* allcolumns */ 1130 lens_i[row]++; /* lens_i[row] += max2 ? */ 1131 } 1132 } 1133 } 1134 } 1135 } 1136 } 1137 ierr = PetscFree(r_waits3);CHKERRQ(ierr); 1138 if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits3,s_status3);CHKERRQ(ierr);} 1139 ierr = PetscFree2(r_status3,s_status3);CHKERRQ(ierr); 1140 ierr = PetscFree(s_waits3);CHKERRQ(ierr); 1141 1142 /* Create the submatrices */ 1143 for (i=0; i<ismax; i++) { 1144 PetscInt bs_tmp; 1145 if (ijonly) bs_tmp = 1; 1146 else bs_tmp = bs; 1147 1148 ierr = MatCreate(PETSC_COMM_SELF,submats+i);CHKERRQ(ierr); 1149 ierr = MatSetSizes(submats[i],nrow[i]*bs_tmp,ncol[i]*bs_tmp,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 1150 1151 ierr = MatSetType(submats[i],((PetscObject)A)->type_name);CHKERRQ(ierr); 1152 ierr = MatSeqBAIJSetPreallocation(submats[i],bs_tmp,0,lens[i]);CHKERRQ(ierr); 1153 ierr = MatSeqSBAIJSetPreallocation(submats[i],bs_tmp,0,lens[i]);CHKERRQ(ierr); /* this subroutine is used by SBAIJ routines */ 1154 1155 /* create struct Mat_SubMat and attached it to submat */ 1156 ierr = PetscNew(&smat_i);CHKERRQ(ierr); 1157 subc = (Mat_SeqBAIJ*)submats[i]->data; 1158 subc->submatis1 = smat_i; 1159 smats[i] = smat_i; 1160 1161 smat_i->destroy = submats[i]->ops->destroy; 1162 submats[i]->ops->destroy = MatDestroy_MPIBAIJ_MatGetSubmatrices; 1163 submats[i]->factortype = C->factortype; 1164 1165 smat_i->id = i; 1166 smat_i->nrqs = nrqs; 1167 smat_i->nrqr = nrqr; 1168 smat_i->rbuf1 = rbuf1; 1169 smat_i->rbuf2 = rbuf2; 1170 smat_i->rbuf3 = rbuf3; 1171 smat_i->sbuf2 = sbuf2; 1172 smat_i->req_source2 = req_source2; 1173 1174 smat_i->sbuf1 = sbuf1; 1175 smat_i->ptr = ptr; 1176 smat_i->tmp = tmp; 1177 smat_i->ctr = ctr; 1178 1179 smat_i->pa = pa; 1180 smat_i->req_size = req_size; 1181 smat_i->req_source1 = req_source1; 1182 1183 smat_i->allcolumns = allcolumns[i]; 1184 smat_i->singleis = PETSC_FALSE; 1185 smat_i->row2proc = row2proc[i]; 1186 smat_i->rmap = rmap[i]; 1187 smat_i->cmap = cmap[i]; 1188 } 1189 1190 if (ismax) {ierr = PetscFree(lens[0]);CHKERRQ(ierr);} 1191 ierr = PetscFree(lens);CHKERRQ(ierr); 1192 ierr = PetscFree(sbuf_aj[0]);CHKERRQ(ierr); 1193 ierr = PetscFree(sbuf_aj);CHKERRQ(ierr); 1194 1195 } /* endof scall == MAT_INITIAL_MATRIX */ 1196 1197 /* Post recv matrix values */ 1198 ierr = PetscObjectGetNewTag((PetscObject)C,&tag4);CHKERRQ(ierr); 1199 ierr = PetscMalloc1(nrqs+1,&rbuf4);CHKERRQ(ierr); 1200 ierr = PetscMalloc1(nrqs+1,&r_waits4);CHKERRQ(ierr); 1201 ierr = PetscMalloc1(nrqs+1,&r_status4);CHKERRQ(ierr); 1202 ierr = PetscMalloc1(nrqr+1,&s_status4);CHKERRQ(ierr); 1203 for (i=0; i<nrqs; ++i) { 1204 ierr = PetscMalloc1(rbuf2[i][0]*bs2,&rbuf4[i]);CHKERRQ(ierr); 1205 ierr = MPI_Irecv(rbuf4[i],rbuf2[i][0]*bs2,MPIU_SCALAR,req_source2[i],tag4,comm,r_waits4+i);CHKERRQ(ierr); 1206 } 1207 1208 /* Allocate sending buffers for a->a, and send them off */ 1209 ierr = PetscMalloc1(nrqr+1,&sbuf_aa);CHKERRQ(ierr); 1210 for (i=0,j=0; i<nrqr; i++) j += req_size[i]; 1211 1212 ierr = PetscMalloc1((j+1)*bs2,&sbuf_aa[0]);CHKERRQ(ierr); 1213 for (i=1; i<nrqr; i++) sbuf_aa[i] = sbuf_aa[i-1] + req_size[i-1]*bs2; 1214 1215 ierr = PetscMalloc1(nrqr+1,&s_waits4);CHKERRQ(ierr); 1216 1217 for (i=0; i<nrqr; i++) { 1218 rbuf1_i = rbuf1[i]; 1219 sbuf_aa_i = sbuf_aa[i]; 1220 ct1 = 2*rbuf1_i[0]+1; 1221 ct2 = 0; 1222 for (j=1,max1=rbuf1_i[0]; j<=max1; j++) { 1223 kmax = rbuf1_i[2*j]; 1224 for (k=0; k<kmax; k++,ct1++) { 1225 row = rbuf1_i[ct1] - rstart; 1226 nzA = a_i[row+1] - a_i[row]; 1227 nzB = b_i[row+1] - b_i[row]; 1228 ncols = nzA + nzB; 1229 cworkB = b_j + b_i[row]; 1230 vworkA = a_a + a_i[row]*bs2; 1231 vworkB = b_a + b_i[row]*bs2; 1232 1233 /* load the column values for this row into vals*/ 1234 vals = sbuf_aa_i+ct2*bs2; 1235 for (l=0; l<nzB; l++) { 1236 if ((bmap[cworkB[l]]) < cstart) { 1237 ierr = PetscMemcpy(vals+l*bs2,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1238 } else break; 1239 } 1240 imark = l; 1241 for (l=0; l<nzA; l++) { 1242 ierr = PetscMemcpy(vals+(imark+l)*bs2,vworkA+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);//error in proc[1]??? 1243 } 1244 for (l=imark; l<nzB; l++) { 1245 ierr = PetscMemcpy(vals+(nzA+l)*bs2,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1246 } 1247 1248 ct2 += ncols; 1249 } 1250 } 1251 ierr = MPI_Isend(sbuf_aa_i,req_size[i]*bs2,MPIU_SCALAR,req_source1[i],tag4,comm,s_waits4+i);CHKERRQ(ierr); 1252 } 1253 1254 if (!ismax) { 1255 ierr = PetscFree(rbuf1[0]);CHKERRQ(ierr); 1256 ierr = PetscFree(rbuf1);CHKERRQ(ierr); 1257 } 1258 1259 /* Assemble the matrices */ 1260 /* First assemble the local rows */ 1261 for (i=0; i<ismax; i++) { 1262 row2proc_i = row2proc[i]; 1263 subc = (Mat_SeqBAIJ*)submats[i]->data; 1264 imat_ilen = subc->ilen; 1265 imat_j = subc->j; 1266 imat_i = subc->i; 1267 imat_a = subc->a; 1268 1269 if (!allcolumns[i]) cmap_i = cmap[i]; 1270 rmap_i = rmap[i]; 1271 irow_i = irow[i]; 1272 jmax = nrow[i]; 1273 for (j=0; j<jmax; j++) { 1274 row = irow_i[j]; 1275 proc = row2proc_i[j]; 1276 1277 if (proc == rank) { 1278 1279 row = row - rstart; 1280 nzA = a_i[row+1] - a_i[row]; 1281 nzB = b_i[row+1] - b_i[row]; 1282 cworkA = a_j + a_i[row]; 1283 cworkB = b_j + b_i[row]; 1284 if (!ijonly) { 1285 vworkA = a_a + a_i[row]*bs2; 1286 vworkB = b_a + b_i[row]*bs2; 1287 } 1288 #if defined(PETSC_USE_CTABLE) 1289 ierr = PetscTableFind(rmap_i,row+rstart+1,&row);CHKERRQ(ierr); 1290 row--; 1291 if (row < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"row not found in table"); 1292 #else 1293 row = rmap_i[row + rstart]; 1294 #endif 1295 mat_i = imat_i[row]; 1296 if (!ijonly) mat_a = imat_a + mat_i*bs2; 1297 mat_j = imat_j + mat_i; 1298 ilen_row = imat_ilen[row]; 1299 1300 /* load the column indices for this row into cols*/ 1301 if (!allcolumns[i]) { 1302 for (l=0; l<nzB; l++) { 1303 if ((ctmp = bmap[cworkB[l]]) < cstart) { 1304 #if defined(PETSC_USE_CTABLE) 1305 ierr = PetscTableFind(cmap_i,ctmp+1,&tcol);CHKERRQ(ierr); 1306 if (tcol) { 1307 #else 1308 if ((tcol = cmap_i[ctmp])) { 1309 #endif 1310 *mat_j++ = tcol - 1; 1311 ierr = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1312 mat_a += bs2; 1313 ilen_row++; 1314 } 1315 } else break; 1316 } 1317 imark = l; 1318 for (l=0; l<nzA; l++) { 1319 #if defined(PETSC_USE_CTABLE) 1320 ierr = PetscTableFind(cmap_i,cstart+cworkA[l]+1,&tcol);CHKERRQ(ierr); 1321 if (tcol) { 1322 #else 1323 if ((tcol = cmap_i[cstart + cworkA[l]])) { 1324 #endif 1325 *mat_j++ = tcol - 1; 1326 if (!ijonly) { 1327 ierr = PetscMemcpy(mat_a,vworkA+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1328 mat_a += bs2; 1329 } 1330 ilen_row++; 1331 } 1332 } 1333 for (l=imark; l<nzB; l++) { 1334 #if defined(PETSC_USE_CTABLE) 1335 ierr = PetscTableFind(cmap_i,bmap[cworkB[l]]+1,&tcol);CHKERRQ(ierr); 1336 if (tcol) { 1337 #else 1338 if ((tcol = cmap_i[bmap[cworkB[l]]])) { 1339 #endif 1340 *mat_j++ = tcol - 1; 1341 if (!ijonly) { 1342 ierr = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1343 mat_a += bs2; 1344 } 1345 ilen_row++; 1346 } 1347 } 1348 } else { /* allcolumns */ 1349 for (l=0; l<nzB; l++) { 1350 if ((ctmp = bmap[cworkB[l]]) < cstart) { 1351 *mat_j++ = ctmp; 1352 ierr = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1353 mat_a += bs2; 1354 ilen_row++; 1355 } else break; 1356 } 1357 imark = l; 1358 for (l=0; l<nzA; l++) { 1359 *mat_j++ = cstart+cworkA[l]; 1360 if (!ijonly) { 1361 ierr = PetscMemcpy(mat_a,vworkA+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1362 mat_a += bs2; 1363 } 1364 ilen_row++; 1365 } 1366 for (l=imark; l<nzB; l++) { 1367 *mat_j++ = bmap[cworkB[l]]; 1368 if (!ijonly) { 1369 ierr = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1370 mat_a += bs2; 1371 } 1372 ilen_row++; 1373 } 1374 } 1375 imat_ilen[row] = ilen_row; 1376 } 1377 } 1378 } 1379 1380 /* Now assemble the off proc rows */ 1381 if (!ijonly) { 1382 ierr = MPI_Waitall(nrqs,r_waits4,r_status4);CHKERRQ(ierr); 1383 } 1384 for (tmp2=0; tmp2<nrqs; tmp2++) { 1385 sbuf1_i = sbuf1[pa[tmp2]]; 1386 jmax = sbuf1_i[0]; 1387 ct1 = 2*jmax + 1; 1388 ct2 = 0; 1389 rbuf2_i = rbuf2[tmp2]; 1390 rbuf3_i = rbuf3[tmp2]; 1391 if (!ijonly) rbuf4_i = rbuf4[tmp2]; 1392 for (j=1; j<=jmax; j++) { 1393 is_no = sbuf1_i[2*j-1]; 1394 rmap_i = rmap[is_no]; 1395 if (!allcolumns[is_no]) cmap_i = cmap[is_no]; 1396 subc = (Mat_SeqBAIJ*)submats[is_no]->data; 1397 imat_ilen = subc->ilen; 1398 imat_j = subc->j; 1399 imat_i = subc->i; 1400 if (!ijonly) imat_a = subc->a; 1401 max1 = sbuf1_i[2*j]; 1402 for (k=0; k<max1; k++,ct1++) { /* for each recved block row */ 1403 row = sbuf1_i[ct1]; 1404 #if defined(PETSC_USE_CTABLE) 1405 ierr = PetscTableFind(rmap_i,row+1,&row);CHKERRQ(ierr); 1406 row--; 1407 if (row < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"row not found in table"); 1408 #else 1409 row = rmap_i[row]; 1410 #endif 1411 ilen = imat_ilen[row]; 1412 mat_i = imat_i[row]; 1413 if (!ijonly) mat_a = imat_a + mat_i*bs2; 1414 mat_j = imat_j + mat_i; 1415 max2 = rbuf2_i[ct1]; 1416 if (!allcolumns[is_no]) { 1417 for (l=0; l<max2; l++,ct2++) { 1418 #if defined(PETSC_USE_CTABLE) 1419 ierr = PetscTableFind(cmap_i,rbuf3_i[ct2]+1,&tcol);CHKERRQ(ierr); 1420 #else 1421 tcol = cmap_i[rbuf3_i[ct2]]; 1422 #endif 1423 if (tcol) { 1424 *mat_j++ = tcol - 1; 1425 if (!ijonly) { 1426 ierr = PetscMemcpy(mat_a,rbuf4_i+ct2*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1427 mat_a += bs2; 1428 } 1429 //*mat_a++ = rbuf4_i[ct2]; 1430 ilen++; 1431 } 1432 } 1433 } else { /* allcolumns */ 1434 for (l=0; l<max2; l++,ct2++) { 1435 *mat_j++ = rbuf3_i[ct2]; /* same global column index of C */ 1436 //*mat_a++ = rbuf4_i[ct2]; 1437 if (!ijonly) { 1438 ierr = PetscMemcpy(mat_a,rbuf4_i+ct2*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1439 mat_a += bs2; 1440 } 1441 ilen++; 1442 } 1443 } 1444 imat_ilen[row] = ilen; 1445 } 1446 } 1447 } 1448 1449 if (!iscsorted) { /* sort column indices of the rows -- not done yet!!! */ 1450 for (i=0; i<ismax; i++) { 1451 subc = (Mat_SeqBAIJ*)submats[i]->data; 1452 imat_j = subc->j; 1453 imat_i = subc->i; 1454 imat_a = subc->a; 1455 imat_ilen = subc->ilen; 1456 1457 if (allcolumns[i]) continue; 1458 jmax = nrow[i]; 1459 for (j=0; j<jmax; j++) { 1460 PetscInt ilen; 1461 1462 mat_i = imat_i[j]; 1463 if (!ijonly) mat_a = imat_a + mat_i; 1464 mat_j = imat_j + mat_i; 1465 ilen = imat_ilen[j]; 1466 ierr = PetscSortIntWithScalarArray(ilen,mat_j,mat_a);CHKERRQ(ierr); 1467 } 1468 } 1469 } 1470 1471 ierr = PetscFree(r_status4);CHKERRQ(ierr); 1472 ierr = PetscFree(r_waits4);CHKERRQ(ierr); 1473 if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits4,s_status4);CHKERRQ(ierr);} 1474 ierr = PetscFree(s_waits4);CHKERRQ(ierr); 1475 ierr = PetscFree(s_status4);CHKERRQ(ierr); 1476 1477 /* Restore the indices */ 1478 for (i=0; i<ismax; i++) { 1479 ierr = ISRestoreIndices(isrow[i],irow+i);CHKERRQ(ierr); 1480 if (!allcolumns[i]) { 1481 ierr = ISRestoreIndices(iscol[i],icol+i);CHKERRQ(ierr); 1482 } 1483 } 1484 1485 for (i=0; i<ismax; i++) { 1486 ierr = MatAssemblyBegin(submats[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1487 ierr = MatAssemblyEnd(submats[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1488 } 1489 1490 /* Destroy allocated memory */ 1491 if (!ismax) { 1492 ierr = PetscFree(pa);CHKERRQ(ierr); 1493 1494 ierr = PetscFree4(sbuf1,ptr,tmp,ctr);CHKERRQ(ierr); 1495 for (i=0; i<nrqr; ++i) { 1496 ierr = PetscFree(sbuf2[i]);CHKERRQ(ierr); 1497 } 1498 for (i=0; i<nrqs; ++i) { 1499 ierr = PetscFree(rbuf3[i]);CHKERRQ(ierr); 1500 } 1501 1502 ierr = PetscFree3(sbuf2,req_size,req_source1);CHKERRQ(ierr); 1503 ierr = PetscFree3(req_source2,rbuf2,rbuf3);CHKERRQ(ierr); 1504 } 1505 1506 ierr = PetscFree(sbuf_aa[0]);CHKERRQ(ierr); 1507 ierr = PetscFree(sbuf_aa);CHKERRQ(ierr); 1508 ierr = PetscFree5(irow,icol,nrow,ncol,issorted);CHKERRQ(ierr); 1509 1510 for (i=0; i<nrqs; ++i) { 1511 ierr = PetscFree(rbuf4[i]);CHKERRQ(ierr); 1512 } 1513 ierr = PetscFree(rbuf4);CHKERRQ(ierr); 1514 1515 ierr = PetscFree5(smats,row2proc,cmap,rmap,allcolumns);CHKERRQ(ierr); 1516 PetscFunctionReturn(0); 1517 } 1518 1519 //------------ endof new ------------- 1520 1521 PetscErrorCode MatGetSubMatrices_MPIBAIJ_local(Mat C,PetscInt ismax,const IS isrow[],const IS iscol[],MatReuse scall,PetscBool *allrows,PetscBool *allcolumns,Mat *submats) 1522 { 1523 Mat_MPIBAIJ *c = (Mat_MPIBAIJ*)C->data; 1524 Mat A = c->A; 1525 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)c->B->data,*mat; 1526 const PetscInt **irow,**icol,*irow_i; 1527 PetscInt *nrow,*ncol,*w3,*w4,start; 1528 PetscErrorCode ierr; 1529 PetscMPIInt size,tag0,tag1,tag2,tag3,*w1,*w2,nrqr,idex,end,proc; 1530 PetscInt **sbuf1,**sbuf2,rank,i,j,k,l,ct1,ct2,**rbuf1,row; 1531 PetscInt nrqs,msz,**ptr,*req_size,*ctr,*pa,*tmp,tcol; 1532 PetscInt **rbuf3,*req_source,**sbuf_aj,**rbuf2,max1,max2; 1533 PetscInt **lens,is_no,ncols,*cols,mat_i,*mat_j,tmp2,jmax; 1534 PetscInt ctr_j,*sbuf1_j,*sbuf_aj_i,*rbuf1_i,kmax,*lens_i; 1535 PetscInt bs =C->rmap->bs,bs2=c->bs2,*a_j=a->j,*b_j=b->j,*cworkA,*cworkB; 1536 PetscInt cstart = c->cstartbs,nzA,nzB,*a_i=a->i,*b_i=b->i,imark; 1537 PetscInt *bmap = c->garray,ctmp,rstart=c->rstartbs; 1538 MPI_Request *s_waits1,*r_waits1,*s_waits2,*r_waits2,*r_waits3,*s_waits3; 1539 MPI_Status *r_status1,*r_status2,*s_status1,*s_status3,*s_status2,*r_status3; 1540 MPI_Comm comm; 1541 PetscBool flag; 1542 PetscMPIInt *onodes1,*olengths1; 1543 PetscBool ijonly=c->ijonly; /* private flag indicates only matrix data structures are requested */ 1544 1545 /* variables below are used for the matrix numerical values - case of !ijonly */ 1546 MPI_Request *r_waits4,*s_waits4; 1547 MPI_Status *r_status4,*s_status4; 1548 MatScalar **rbuf4,**sbuf_aa,*vals,*mat_a = NULL,*sbuf_aa_i,*vworkA = NULL,*vworkB = NULL; 1549 MatScalar *a_a=a->a,*b_a=b->a; 1550 1551 #if defined(PETSC_USE_CTABLE) 1552 PetscInt tt; 1553 PetscTable *rmap,*cmap,rmap_i,cmap_i=NULL; 1554 #else 1555 PetscInt **cmap,*cmap_i=NULL,*rtable,*rmap_i,**rmap, Mbs = c->Mbs; 1556 #endif 1557 1558 PetscFunctionBegin; 1559 ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr); 1560 tag0 = ((PetscObject)C)->tag; 1561 size = c->size; 1562 rank = c->rank; 1563 //if (!rank) printf("MatGetSubMatrices_MPIBAIJ scall %d, bs %d\n",scall,bs); 1564 1565 /* Get some new tags to keep the communication clean */ 1566 ierr = PetscObjectGetNewTag((PetscObject)C,&tag1);CHKERRQ(ierr); 1567 ierr = PetscObjectGetNewTag((PetscObject)C,&tag2);CHKERRQ(ierr); 1568 ierr = PetscObjectGetNewTag((PetscObject)C,&tag3);CHKERRQ(ierr); 1569 1570 #if defined(PETSC_USE_CTABLE) 1571 ierr = PetscMalloc4(ismax,&irow,ismax,&icol,ismax,&nrow,ismax,&ncol);CHKERRQ(ierr); 1572 #else 1573 ierr = PetscMalloc5(ismax,&irow,ismax,&icol,ismax,&nrow,ismax,&ncol,Mbs+1,&rtable);CHKERRQ(ierr); 1574 /* Create hash table for the mapping :row -> proc*/ 1575 for (i=0,j=0; i<size; i++) { 1576 jmax = C->rmap->range[i+1]/bs; 1577 for (; j<jmax; j++) rtable[j] = i; 1578 } 1579 #endif 1580 1581 for (i=0; i<ismax; i++) { 1582 if (allrows[i]) { 1583 irow[i] = NULL; 1584 nrow[i] = C->rmap->N/bs; 1585 } else { 1586 ierr = ISGetIndices(isrow[i],&irow[i]);CHKERRQ(ierr); 1587 ierr = ISGetLocalSize(isrow[i],&nrow[i]);CHKERRQ(ierr); 1588 } 1589 1590 if (allcolumns[i]) { 1591 icol[i] = NULL; 1592 ncol[i] = C->cmap->N/bs; 1593 } else { 1594 ierr = ISGetIndices(iscol[i],&icol[i]);CHKERRQ(ierr); 1595 ierr = ISGetLocalSize(iscol[i],&ncol[i]);CHKERRQ(ierr); 1596 } 1597 } 1598 1599 /* evaluate communication - mesg to who,length of mesg,and buffer space 1600 required. Based on this, buffers are allocated, and data copied into them*/ 1601 ierr = PetscCalloc4(size,&w1,size,&w2,size,&w3,size,&w4);CHKERRQ(ierr); 1602 for (i=0; i<ismax; i++) { 1603 ierr = PetscMemzero(w4,size*sizeof(PetscInt));CHKERRQ(ierr); /* initialise work vector*/ 1604 jmax = nrow[i]; 1605 irow_i = irow[i]; 1606 for (j=0; j<jmax; j++) { 1607 if (allrows[i]) row = j; 1608 else row = irow_i[j]; 1609 1610 #if defined(PETSC_USE_CTABLE) 1611 ierr = PetscGetProc(row,size,c->rangebs,&proc);CHKERRQ(ierr); 1612 #else 1613 proc = rtable[row]; 1614 #endif 1615 w4[proc]++; 1616 } 1617 for (j=0; j<size; j++) { 1618 if (w4[j]) { w1[j] += w4[j]; w3[j]++;} 1619 } 1620 } 1621 1622 nrqs = 0; /* no of outgoing messages */ 1623 msz = 0; /* total mesg length for all proc */ 1624 w1[rank] = 0; /* no mesg sent to intself */ 1625 w3[rank] = 0; 1626 for (i=0; i<size; i++) { 1627 if (w1[i]) { w2[i] = 1; nrqs++;} /* there exists a message to proc i */ 1628 } 1629 ierr = PetscMalloc1(nrqs+1,&pa);CHKERRQ(ierr); /*(proc -array)*/ 1630 for (i=0,j=0; i<size; i++) { 1631 if (w1[i]) { pa[j] = i; j++; } 1632 } 1633 1634 /* Each message would have a header = 1 + 2*(no of IS) + data */ 1635 for (i=0; i<nrqs; i++) { 1636 j = pa[i]; 1637 w1[j] += w2[j] + 2* w3[j]; 1638 msz += w1[j]; 1639 } 1640 1641 /* Determine the number of messages to expect, their lengths, from from-ids */ 1642 ierr = PetscGatherNumberOfMessages(comm,w2,w1,&nrqr);CHKERRQ(ierr); 1643 ierr = PetscGatherMessageLengths(comm,nrqs,nrqr,w1,&onodes1,&olengths1);CHKERRQ(ierr); 1644 1645 /* Now post the Irecvs corresponding to these messages */ 1646 ierr = PetscPostIrecvInt(comm,tag0,nrqr,onodes1,olengths1,&rbuf1,&r_waits1);CHKERRQ(ierr); 1647 1648 ierr = PetscFree(onodes1);CHKERRQ(ierr); 1649 ierr = PetscFree(olengths1);CHKERRQ(ierr); 1650 1651 /* Allocate Memory for outgoing messages */ 1652 ierr = PetscMalloc4(size,&sbuf1,size,&ptr,2*msz,&tmp,size,&ctr);CHKERRQ(ierr); 1653 ierr = PetscMemzero(sbuf1,size*sizeof(PetscInt*));CHKERRQ(ierr); 1654 ierr = PetscMemzero(ptr,size*sizeof(PetscInt*));CHKERRQ(ierr); 1655 { 1656 PetscInt *iptr = tmp,ict = 0; 1657 for (i=0; i<nrqs; i++) { 1658 j = pa[i]; 1659 iptr += ict; 1660 sbuf1[j] = iptr; 1661 ict = w1[j]; 1662 } 1663 } 1664 1665 /* Form the outgoing messages */ 1666 /* Initialise the header space */ 1667 for (i=0; i<nrqs; i++) { 1668 j = pa[i]; 1669 sbuf1[j][0] = 0; 1670 ierr = PetscMemzero(sbuf1[j]+1,2*w3[j]*sizeof(PetscInt));CHKERRQ(ierr); 1671 ptr[j] = sbuf1[j] + 2*w3[j] + 1; 1672 } 1673 1674 /* Parse the isrow and copy data into outbuf */ 1675 for (i=0; i<ismax; i++) { 1676 ierr = PetscMemzero(ctr,size*sizeof(PetscInt));CHKERRQ(ierr); 1677 irow_i = irow[i]; 1678 jmax = nrow[i]; 1679 for (j=0; j<jmax; j++) { /* parse the indices of each IS */ 1680 if (allrows[i]) row = j; 1681 else row = irow_i[j]; 1682 1683 #if defined(PETSC_USE_CTABLE) 1684 ierr = PetscGetProc(row,size,c->rangebs,&proc);CHKERRQ(ierr); 1685 #else 1686 proc = rtable[row]; 1687 #endif 1688 if (proc != rank) { /* copy to the outgoing buf*/ 1689 ctr[proc]++; 1690 *ptr[proc] = row; 1691 ptr[proc]++; 1692 } 1693 } 1694 /* Update the headers for the current IS */ 1695 for (j=0; j<size; j++) { /* Can Optimise this loop too */ 1696 if ((ctr_j = ctr[j])) { 1697 sbuf1_j = sbuf1[j]; 1698 k = ++sbuf1_j[0]; 1699 sbuf1_j[2*k] = ctr_j; 1700 sbuf1_j[2*k-1] = i; 1701 } 1702 } 1703 } 1704 1705 /* Now post the sends */ 1706 ierr = PetscMalloc1(nrqs+1,&s_waits1);CHKERRQ(ierr); 1707 for (i=0; i<nrqs; ++i) { 1708 j = pa[i]; 1709 ierr = MPI_Isend(sbuf1[j],w1[j],MPIU_INT,j,tag0,comm,s_waits1+i);CHKERRQ(ierr); 1710 } 1711 1712 /* Post Recieves to capture the buffer size */ 1713 ierr = PetscMalloc1(nrqs+1,&r_waits2);CHKERRQ(ierr); 1714 ierr = PetscMalloc1(nrqs+1,&rbuf2);CHKERRQ(ierr); 1715 rbuf2[0] = tmp + msz; 1716 for (i=1; i<nrqs; ++i) { 1717 rbuf2[i] = rbuf2[i-1]+w1[pa[i-1]]; 1718 } 1719 for (i=0; i<nrqs; ++i) { 1720 j = pa[i]; 1721 ierr = MPI_Irecv(rbuf2[i],w1[j],MPIU_INT,j,tag1,comm,r_waits2+i);CHKERRQ(ierr); 1722 } 1723 1724 /* Send to other procs the buf size they should allocate */ 1725 1726 /* Receive messages*/ 1727 ierr = PetscMalloc1(nrqr+1,&s_waits2);CHKERRQ(ierr); 1728 ierr = PetscMalloc1(nrqr+1,&r_status1);CHKERRQ(ierr); 1729 ierr = PetscMalloc3(nrqr+1,&sbuf2,nrqr,&req_size,nrqr,&req_source);CHKERRQ(ierr); 1730 { 1731 Mat_SeqBAIJ *sA = (Mat_SeqBAIJ*)c->A->data,*sB = (Mat_SeqBAIJ*)c->B->data; 1732 PetscInt *sAi = sA->i,*sBi = sB->i,id,*sbuf2_i; 1733 1734 for (i=0; i<nrqr; ++i) { 1735 ierr = MPI_Waitany(nrqr,r_waits1,&idex,r_status1+i);CHKERRQ(ierr); 1736 1737 req_size[idex] = 0; 1738 rbuf1_i = rbuf1[idex]; 1739 start = 2*rbuf1_i[0] + 1; 1740 ierr = MPI_Get_count(r_status1+i,MPIU_INT,&end);CHKERRQ(ierr); 1741 ierr = PetscMalloc1(end,&sbuf2[idex]);CHKERRQ(ierr); 1742 sbuf2_i = sbuf2[idex]; 1743 for (j=start; j<end; j++) { 1744 id = rbuf1_i[j] - rstart; 1745 ncols = sAi[id+1] - sAi[id] + sBi[id+1] - sBi[id]; 1746 sbuf2_i[j] = ncols; 1747 req_size[idex] += ncols; 1748 } 1749 req_source[idex] = r_status1[i].MPI_SOURCE; 1750 /* form the header */ 1751 sbuf2_i[0] = req_size[idex]; 1752 for (j=1; j<start; j++) sbuf2_i[j] = rbuf1_i[j]; 1753 ierr = MPI_Isend(sbuf2_i,end,MPIU_INT,req_source[idex],tag1,comm,s_waits2+i);CHKERRQ(ierr); 1754 } 1755 } 1756 ierr = PetscFree(r_status1);CHKERRQ(ierr); 1757 ierr = PetscFree(r_waits1);CHKERRQ(ierr); 1758 1759 /* recv buffer sizes */ 1760 /* Receive messages*/ 1761 ierr = PetscMalloc1(nrqs+1,&rbuf3);CHKERRQ(ierr); 1762 ierr = PetscMalloc1(nrqs+1,&r_waits3);CHKERRQ(ierr); 1763 ierr = PetscMalloc1(nrqs+1,&r_status2);CHKERRQ(ierr); 1764 if (!ijonly) { 1765 ierr = PetscMalloc1(nrqs+1,&rbuf4);CHKERRQ(ierr); 1766 ierr = PetscMalloc1(nrqs+1,&r_waits4);CHKERRQ(ierr); 1767 } 1768 1769 for (i=0; i<nrqs; ++i) { 1770 ierr = MPI_Waitany(nrqs,r_waits2,&idex,r_status2+i);CHKERRQ(ierr); 1771 ierr = PetscMalloc1(rbuf2[idex][0],&rbuf3[idex]);CHKERRQ(ierr); 1772 ierr = MPI_Irecv(rbuf3[idex],rbuf2[idex][0],MPIU_INT,r_status2[i].MPI_SOURCE,tag2,comm,r_waits3+idex);CHKERRQ(ierr); 1773 if (!ijonly) { 1774 ierr = PetscMalloc1(rbuf2[idex][0]*bs2,&rbuf4[idex]);CHKERRQ(ierr); 1775 ierr = MPI_Irecv(rbuf4[idex],rbuf2[idex][0]*bs2,MPIU_MATSCALAR,r_status2[i].MPI_SOURCE,tag3,comm,r_waits4+idex);CHKERRQ(ierr); 1776 } 1777 } 1778 ierr = PetscFree(r_status2);CHKERRQ(ierr); 1779 ierr = PetscFree(r_waits2);CHKERRQ(ierr); 1780 1781 /* Wait on sends1 and sends2 */ 1782 ierr = PetscMalloc1(nrqs+1,&s_status1);CHKERRQ(ierr); 1783 ierr = PetscMalloc1(nrqr+1,&s_status2);CHKERRQ(ierr); 1784 1785 if (nrqs) {ierr = MPI_Waitall(nrqs,s_waits1,s_status1);CHKERRQ(ierr);} 1786 if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits2,s_status2);CHKERRQ(ierr);} 1787 ierr = PetscFree(s_status1);CHKERRQ(ierr); 1788 ierr = PetscFree(s_status2);CHKERRQ(ierr); 1789 ierr = PetscFree(s_waits1);CHKERRQ(ierr); 1790 ierr = PetscFree(s_waits2);CHKERRQ(ierr); 1791 1792 /* Now allocate buffers for a->j, and send them off */ 1793 ierr = PetscMalloc1(nrqr+1,&sbuf_aj);CHKERRQ(ierr); 1794 for (i=0,j=0; i<nrqr; i++) j += req_size[i]; 1795 ierr = PetscMalloc1(j+1,&sbuf_aj[0]);CHKERRQ(ierr); 1796 for (i=1; i<nrqr; i++) sbuf_aj[i] = sbuf_aj[i-1] + req_size[i-1]; 1797 1798 ierr = PetscMalloc1(nrqr+1,&s_waits3);CHKERRQ(ierr); 1799 { 1800 for (i=0; i<nrqr; i++) { 1801 rbuf1_i = rbuf1[i]; 1802 sbuf_aj_i = sbuf_aj[i]; 1803 ct1 = 2*rbuf1_i[0] + 1; 1804 ct2 = 0; 1805 for (j=1,max1=rbuf1_i[0]; j<=max1; j++) { 1806 kmax = rbuf1[i][2*j]; 1807 for (k=0; k<kmax; k++,ct1++) { 1808 row = rbuf1_i[ct1] - rstart; 1809 nzA = a_i[row+1] - a_i[row]; 1810 nzB = b_i[row+1] - b_i[row]; 1811 ncols = nzA + nzB; 1812 cworkA = a_j + a_i[row]; 1813 cworkB = b_j + b_i[row]; 1814 1815 /* load the column indices for this row into cols*/ 1816 cols = sbuf_aj_i + ct2; 1817 for (l=0; l<nzB; l++) { 1818 if ((ctmp = bmap[cworkB[l]]) < cstart) cols[l] = ctmp; 1819 else break; 1820 } 1821 imark = l; 1822 for (l=0; l<nzA; l++) cols[imark+l] = cstart + cworkA[l]; 1823 for (l=imark; l<nzB; l++) cols[nzA+l] = bmap[cworkB[l]]; 1824 ct2 += ncols; 1825 } 1826 } 1827 ierr = MPI_Isend(sbuf_aj_i,req_size[i],MPIU_INT,req_source[i],tag2,comm,s_waits3+i);CHKERRQ(ierr); 1828 } 1829 } 1830 ierr = PetscMalloc1(nrqs+1,&r_status3);CHKERRQ(ierr); 1831 ierr = PetscMalloc1(nrqr+1,&s_status3);CHKERRQ(ierr); 1832 1833 /* Allocate buffers for a->a, and send them off */ 1834 if (!ijonly) { 1835 ierr = PetscMalloc1(nrqr+1,&sbuf_aa);CHKERRQ(ierr); 1836 for (i=0,j=0; i<nrqr; i++) j += req_size[i]; 1837 ierr = PetscMalloc1((j+1)*bs2,&sbuf_aa[0]);CHKERRQ(ierr); 1838 for (i=1; i<nrqr; i++) sbuf_aa[i] = sbuf_aa[i-1] + req_size[i-1]*bs2; 1839 1840 ierr = PetscMalloc1(nrqr+1,&s_waits4);CHKERRQ(ierr); 1841 { 1842 for (i=0; i<nrqr; i++) { 1843 rbuf1_i = rbuf1[i]; 1844 sbuf_aa_i = sbuf_aa[i]; 1845 ct1 = 2*rbuf1_i[0]+1; 1846 ct2 = 0; 1847 for (j=1,max1=rbuf1_i[0]; j<=max1; j++) { 1848 kmax = rbuf1_i[2*j]; 1849 for (k=0; k<kmax; k++,ct1++) { 1850 row = rbuf1_i[ct1] - rstart; 1851 nzA = a_i[row+1] - a_i[row]; 1852 nzB = b_i[row+1] - b_i[row]; 1853 ncols = nzA + nzB; 1854 cworkB = b_j + b_i[row]; 1855 vworkA = a_a + a_i[row]*bs2; 1856 vworkB = b_a + b_i[row]*bs2; 1857 1858 /* load the column values for this row into vals*/ 1859 vals = sbuf_aa_i+ct2*bs2; 1860 for (l=0; l<nzB; l++) { 1861 if ((bmap[cworkB[l]]) < cstart) { 1862 ierr = PetscMemcpy(vals+l*bs2,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1863 } else break; 1864 } 1865 imark = l; 1866 for (l=0; l<nzA; l++) { 1867 ierr = PetscMemcpy(vals+(imark+l)*bs2,vworkA+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1868 } 1869 for (l=imark; l<nzB; l++) { 1870 ierr = PetscMemcpy(vals+(nzA+l)*bs2,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1871 } 1872 ct2 += ncols; 1873 } 1874 } 1875 ierr = MPI_Isend(sbuf_aa_i,req_size[i]*bs2,MPIU_MATSCALAR,req_source[i],tag3,comm,s_waits4+i);CHKERRQ(ierr); 1876 } 1877 } 1878 ierr = PetscMalloc1(nrqs+1,&r_status4);CHKERRQ(ierr); 1879 ierr = PetscMalloc1(nrqr+1,&s_status4);CHKERRQ(ierr); 1880 } 1881 ierr = PetscFree(rbuf1[0]);CHKERRQ(ierr); 1882 ierr = PetscFree(rbuf1);CHKERRQ(ierr); 1883 1884 /* Form the matrix */ 1885 /* create col map: global col of C -> local col of submatrices */ 1886 { 1887 const PetscInt *icol_i; 1888 #if defined(PETSC_USE_CTABLE) 1889 ierr = PetscMalloc1(1+ismax,&cmap);CHKERRQ(ierr); 1890 for (i=0; i<ismax; i++) { 1891 if (!allcolumns[i]) { 1892 ierr = PetscTableCreate(ncol[i]+1,c->Nbs+1,&cmap[i]);CHKERRQ(ierr); 1893 jmax = ncol[i]; 1894 icol_i = icol[i]; 1895 cmap_i = cmap[i]; 1896 for (j=0; j<jmax; j++) { 1897 ierr = PetscTableAdd(cmap_i,icol_i[j]+1,j+1,INSERT_VALUES);CHKERRQ(ierr); 1898 } 1899 } else { 1900 cmap[i] = NULL; 1901 } 1902 } 1903 #else 1904 ierr = PetscMalloc1(ismax,&cmap);CHKERRQ(ierr); 1905 for (i=0; i<ismax; i++) { 1906 if (!allcolumns[i]) { 1907 ierr = PetscCalloc1(c->Nbs,&cmap[i]);CHKERRQ(ierr); 1908 jmax = ncol[i]; 1909 icol_i = icol[i]; 1910 cmap_i = cmap[i]; 1911 for (j=0; j<jmax; j++) { 1912 cmap_i[icol_i[j]] = j+1; 1913 } 1914 } else { /* allcolumns[i] */ 1915 cmap[i] = NULL; 1916 } 1917 } 1918 #endif 1919 } 1920 1921 /* Create lens which is required for MatCreate... */ 1922 for (i=0,j=0; i<ismax; i++) j += nrow[i]; 1923 ierr = PetscMalloc((1+ismax)*sizeof(PetscInt*)+ j*sizeof(PetscInt),&lens);CHKERRQ(ierr); 1924 lens[0] = (PetscInt*)(lens + ismax); 1925 ierr = PetscMemzero(lens[0],j*sizeof(PetscInt));CHKERRQ(ierr); 1926 for (i=1; i<ismax; i++) lens[i] = lens[i-1] + nrow[i-1]; 1927 1928 /* Update lens from local data */ 1929 for (i=0; i<ismax; i++) { 1930 jmax = nrow[i]; 1931 if (!allcolumns[i]) cmap_i = cmap[i]; 1932 irow_i = irow[i]; 1933 lens_i = lens[i]; 1934 for (j=0; j<jmax; j++) { 1935 if (allrows[i]) row = j; 1936 else row = irow_i[j]; 1937 1938 #if defined(PETSC_USE_CTABLE) 1939 ierr = PetscGetProc(row,size,c->rangebs,&proc);CHKERRQ(ierr); 1940 #else 1941 proc = rtable[row]; 1942 #endif 1943 if (proc == rank) { 1944 /* Get indices from matA and then from matB */ 1945 row = row - rstart; 1946 nzA = a_i[row+1] - a_i[row]; 1947 nzB = b_i[row+1] - b_i[row]; 1948 cworkA = a_j + a_i[row]; 1949 cworkB = b_j + b_i[row]; 1950 if (!allcolumns[i]) { 1951 #if defined(PETSC_USE_CTABLE) 1952 for (k=0; k<nzA; k++) { 1953 ierr = PetscTableFind(cmap_i,cstart+cworkA[k]+1,&tt);CHKERRQ(ierr); 1954 if (tt) lens_i[j]++; 1955 } 1956 for (k=0; k<nzB; k++) { 1957 ierr = PetscTableFind(cmap_i,bmap[cworkB[k]]+1,&tt);CHKERRQ(ierr); 1958 if (tt) lens_i[j]++; 1959 } 1960 1961 #else 1962 for (k=0; k<nzA; k++) { 1963 if (cmap_i[cstart + cworkA[k]]) lens_i[j]++; 1964 } 1965 for (k=0; k<nzB; k++) { 1966 if (cmap_i[bmap[cworkB[k]]]) lens_i[j]++; 1967 } 1968 #endif 1969 } else { /* allcolumns */ 1970 lens_i[j] = nzA + nzB; 1971 } 1972 } 1973 } 1974 } 1975 #if defined(PETSC_USE_CTABLE) 1976 /* Create row map*/ 1977 ierr = PetscMalloc1(1+ismax,&rmap);CHKERRQ(ierr); 1978 for (i=0; i<ismax; i++) { 1979 ierr = PetscTableCreate(nrow[i]+1,c->Mbs+1,&rmap[i]);CHKERRQ(ierr); 1980 } 1981 #else 1982 /* Create row map*/ 1983 ierr = PetscMalloc((1+ismax)*sizeof(PetscInt*)+ ismax*Mbs*sizeof(PetscInt),&rmap);CHKERRQ(ierr); 1984 rmap[0] = (PetscInt*)(rmap + ismax); 1985 ierr = PetscMemzero(rmap[0],ismax*Mbs*sizeof(PetscInt));CHKERRQ(ierr); 1986 for (i=1; i<ismax; i++) rmap[i] = rmap[i-1] + Mbs; 1987 #endif 1988 for (i=0; i<ismax; i++) { 1989 irow_i = irow[i]; 1990 jmax = nrow[i]; 1991 #if defined(PETSC_USE_CTABLE) 1992 rmap_i = rmap[i]; 1993 for (j=0; j<jmax; j++) { 1994 if (allrows[i]) { 1995 ierr = PetscTableAdd(rmap_i,j+1,j+1,INSERT_VALUES);CHKERRQ(ierr); 1996 } else { 1997 ierr = PetscTableAdd(rmap_i,irow_i[j]+1,j+1,INSERT_VALUES);CHKERRQ(ierr); 1998 } 1999 } 2000 #else 2001 rmap_i = rmap[i]; 2002 for (j=0; j<jmax; j++) { 2003 if (allrows[i]) rmap_i[j] = j; 2004 else rmap_i[irow_i[j]] = j; 2005 } 2006 #endif 2007 } 2008 2009 /* Update lens from offproc data */ 2010 { 2011 PetscInt *rbuf2_i,*rbuf3_i,*sbuf1_i; 2012 PetscMPIInt ii; 2013 2014 for (tmp2=0; tmp2<nrqs; tmp2++) { 2015 ierr = MPI_Waitany(nrqs,r_waits3,&ii,r_status3+tmp2);CHKERRQ(ierr); 2016 idex = pa[ii]; 2017 sbuf1_i = sbuf1[idex]; 2018 jmax = sbuf1_i[0]; 2019 ct1 = 2*jmax+1; 2020 ct2 = 0; 2021 rbuf2_i = rbuf2[ii]; 2022 rbuf3_i = rbuf3[ii]; 2023 for (j=1; j<=jmax; j++) { 2024 is_no = sbuf1_i[2*j-1]; 2025 max1 = sbuf1_i[2*j]; 2026 lens_i = lens[is_no]; 2027 if (!allcolumns[is_no]) cmap_i = cmap[is_no]; 2028 rmap_i = rmap[is_no]; 2029 for (k=0; k<max1; k++,ct1++) { 2030 #if defined(PETSC_USE_CTABLE) 2031 ierr = PetscTableFind(rmap_i,sbuf1_i[ct1]+1,&row);CHKERRQ(ierr); 2032 row--; 2033 if (row < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"row not found in table"); 2034 #else 2035 row = rmap_i[sbuf1_i[ct1]]; /* the val in the new matrix to be */ 2036 #endif 2037 max2 = rbuf2_i[ct1]; 2038 for (l=0; l<max2; l++,ct2++) { 2039 if (!allcolumns[is_no]) { 2040 #if defined(PETSC_USE_CTABLE) 2041 ierr = PetscTableFind(cmap_i,rbuf3_i[ct2]+1,&tt);CHKERRQ(ierr); 2042 if (tt) lens_i[row]++; 2043 #else 2044 if (cmap_i[rbuf3_i[ct2]]) lens_i[row]++; 2045 #endif 2046 } else { /* allcolumns */ 2047 lens_i[row]++; 2048 } 2049 } 2050 } 2051 } 2052 } 2053 } 2054 ierr = PetscFree(r_status3);CHKERRQ(ierr); 2055 ierr = PetscFree(r_waits3);CHKERRQ(ierr); 2056 if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits3,s_status3);CHKERRQ(ierr);} 2057 ierr = PetscFree(s_status3);CHKERRQ(ierr); 2058 ierr = PetscFree(s_waits3);CHKERRQ(ierr); 2059 2060 /* Create the submatrices */ 2061 if (scall == MAT_REUSE_MATRIX) { 2062 if (ijonly) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP," MAT_REUSE_MATRIX and ijonly is not supported yet"); 2063 /* 2064 Assumes new rows are same length as the old rows, hence bug! 2065 */ 2066 for (i=0; i<ismax; i++) { 2067 mat = (Mat_SeqBAIJ*)(submats[i]->data); 2068 if ((mat->mbs != nrow[i]) || (mat->nbs != ncol[i] || C->rmap->bs != bs)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong size"); 2069 ierr = PetscMemcmp(mat->ilen,lens[i],mat->mbs *sizeof(PetscInt),&flag);CHKERRQ(ierr); 2070 if (!flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Cannot reuse matrix. wrong no of nonzeros"); 2071 /* Initial matrix as if empty */ 2072 ierr = PetscMemzero(mat->ilen,mat->mbs*sizeof(PetscInt));CHKERRQ(ierr); 2073 2074 submats[i]->factortype = C->factortype; 2075 } 2076 } else { 2077 PetscInt bs_tmp; 2078 if (ijonly) bs_tmp = 1; 2079 else bs_tmp = bs; 2080 for (i=0; i<ismax; i++) { 2081 ierr = MatCreate(PETSC_COMM_SELF,submats+i);CHKERRQ(ierr); 2082 ierr = MatSetSizes(submats[i],nrow[i]*bs_tmp,ncol[i]*bs_tmp,nrow[i]*bs_tmp,ncol[i]*bs_tmp);CHKERRQ(ierr); 2083 ierr = MatSetType(submats[i],((PetscObject)A)->type_name);CHKERRQ(ierr); 2084 ierr = MatSeqBAIJSetPreallocation(submats[i],bs_tmp,0,lens[i]);CHKERRQ(ierr); 2085 ierr = MatSeqSBAIJSetPreallocation(submats[i],bs_tmp,0,lens[i]);CHKERRQ(ierr); /* this subroutine is used by SBAIJ routines */ 2086 } 2087 } 2088 2089 /* Assemble the matrices */ 2090 /* First assemble the local rows */ 2091 { 2092 PetscInt ilen_row,*imat_ilen,*imat_j,*imat_i; 2093 MatScalar *imat_a = NULL; 2094 2095 for (i=0; i<ismax; i++) { 2096 mat = (Mat_SeqBAIJ*)submats[i]->data; 2097 imat_ilen = mat->ilen; 2098 imat_j = mat->j; 2099 imat_i = mat->i; 2100 if (!ijonly) imat_a = mat->a; 2101 if (!allcolumns[i]) cmap_i = cmap[i]; 2102 rmap_i = rmap[i]; 2103 irow_i = irow[i]; 2104 jmax = nrow[i]; 2105 for (j=0; j<jmax; j++) { 2106 if (allrows[i]) row = j; 2107 else row = irow_i[j]; 2108 2109 #if defined(PETSC_USE_CTABLE) 2110 ierr = PetscGetProc(row,size,c->rangebs,&proc);CHKERRQ(ierr); 2111 #else 2112 proc = rtable[row]; 2113 #endif 2114 if (proc == rank) { 2115 row = row - rstart; 2116 nzA = a_i[row+1] - a_i[row]; 2117 nzB = b_i[row+1] - b_i[row]; 2118 cworkA = a_j + a_i[row]; 2119 cworkB = b_j + b_i[row]; 2120 if (!ijonly) { 2121 vworkA = a_a + a_i[row]*bs2; 2122 vworkB = b_a + b_i[row]*bs2; 2123 } 2124 #if defined(PETSC_USE_CTABLE) 2125 ierr = PetscTableFind(rmap_i,row+rstart+1,&row);CHKERRQ(ierr); 2126 row--; 2127 if (row < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"row not found in table"); 2128 #else 2129 row = rmap_i[row + rstart]; 2130 #endif 2131 mat_i = imat_i[row]; 2132 if (!ijonly) mat_a = imat_a + mat_i*bs2; 2133 mat_j = imat_j + mat_i; 2134 ilen_row = imat_ilen[row]; 2135 2136 /* load the column indices for this row into cols*/ 2137 if (!allcolumns[i]) { 2138 for (l=0; l<nzB; l++) { 2139 if ((ctmp = bmap[cworkB[l]]) < cstart) { 2140 #if defined(PETSC_USE_CTABLE) 2141 ierr = PetscTableFind(cmap_i,ctmp+1,&tcol);CHKERRQ(ierr); 2142 if (tcol) { 2143 #else 2144 if ((tcol = cmap_i[ctmp])) { 2145 #endif 2146 *mat_j++ = tcol - 1; 2147 ierr = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 2148 mat_a += bs2; 2149 ilen_row++; 2150 } 2151 } else break; 2152 } 2153 imark = l; 2154 for (l=0; l<nzA; l++) { 2155 #if defined(PETSC_USE_CTABLE) 2156 ierr = PetscTableFind(cmap_i,cstart+cworkA[l]+1,&tcol);CHKERRQ(ierr); 2157 if (tcol) { 2158 #else 2159 if ((tcol = cmap_i[cstart + cworkA[l]])) { 2160 #endif 2161 *mat_j++ = tcol - 1; 2162 if (!ijonly) { 2163 ierr = PetscMemcpy(mat_a,vworkA+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 2164 mat_a += bs2; 2165 } 2166 ilen_row++; 2167 } 2168 } 2169 for (l=imark; l<nzB; l++) { 2170 #if defined(PETSC_USE_CTABLE) 2171 ierr = PetscTableFind(cmap_i,bmap[cworkB[l]]+1,&tcol);CHKERRQ(ierr); 2172 if (tcol) { 2173 #else 2174 if ((tcol = cmap_i[bmap[cworkB[l]]])) { 2175 #endif 2176 *mat_j++ = tcol - 1; 2177 if (!ijonly) { 2178 ierr = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 2179 mat_a += bs2; 2180 } 2181 ilen_row++; 2182 } 2183 } 2184 } else { /* allcolumns */ 2185 for (l=0; l<nzB; l++) { 2186 if ((ctmp = bmap[cworkB[l]]) < cstart) { 2187 *mat_j++ = ctmp; 2188 ierr = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 2189 mat_a += bs2; 2190 ilen_row++; 2191 } else break; 2192 } 2193 imark = l; 2194 for (l=0; l<nzA; l++) { 2195 *mat_j++ = cstart+cworkA[l]; 2196 if (!ijonly) { 2197 ierr = PetscMemcpy(mat_a,vworkA+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 2198 mat_a += bs2; 2199 } 2200 ilen_row++; 2201 } 2202 for (l=imark; l<nzB; l++) { 2203 *mat_j++ = bmap[cworkB[l]]; 2204 if (!ijonly) { 2205 ierr = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 2206 mat_a += bs2; 2207 } 2208 ilen_row++; 2209 } 2210 } 2211 imat_ilen[row] = ilen_row; 2212 } 2213 } 2214 } 2215 } 2216 2217 /* Now assemble the off proc rows*/ 2218 { 2219 PetscInt *sbuf1_i,*rbuf2_i,*rbuf3_i,*imat_ilen,ilen; 2220 PetscInt *imat_j,*imat_i; 2221 MatScalar *imat_a = NULL,*rbuf4_i = NULL; 2222 PetscMPIInt ii; 2223 2224 for (tmp2=0; tmp2<nrqs; tmp2++) { 2225 if (ijonly) ii = tmp2; 2226 else { 2227 ierr = MPI_Waitany(nrqs,r_waits4,&ii,r_status4+tmp2);CHKERRQ(ierr); 2228 } 2229 idex = pa[ii]; 2230 sbuf1_i = sbuf1[idex]; 2231 jmax = sbuf1_i[0]; 2232 ct1 = 2*jmax + 1; 2233 ct2 = 0; 2234 rbuf2_i = rbuf2[ii]; 2235 rbuf3_i = rbuf3[ii]; 2236 if (!ijonly) rbuf4_i = rbuf4[ii]; 2237 for (j=1; j<=jmax; j++) { 2238 is_no = sbuf1_i[2*j-1]; 2239 if (!allcolumns[is_no]) cmap_i = cmap[is_no]; 2240 rmap_i = rmap[is_no]; 2241 mat = (Mat_SeqBAIJ*)submats[is_no]->data; 2242 imat_ilen = mat->ilen; 2243 imat_j = mat->j; 2244 imat_i = mat->i; 2245 if (!ijonly) imat_a = mat->a; 2246 max1 = sbuf1_i[2*j]; 2247 for (k=0; k<max1; k++,ct1++) { 2248 row = sbuf1_i[ct1]; 2249 #if defined(PETSC_USE_CTABLE) 2250 ierr = PetscTableFind(rmap_i,row+1,&row);CHKERRQ(ierr); 2251 row--; 2252 if (row < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"row not found in table"); 2253 #else 2254 row = rmap_i[row]; 2255 #endif 2256 ilen = imat_ilen[row]; 2257 mat_i = imat_i[row]; 2258 if (!ijonly) mat_a = imat_a + mat_i*bs2; 2259 mat_j = imat_j + mat_i; 2260 max2 = rbuf2_i[ct1]; 2261 2262 if (!allcolumns[is_no]) { 2263 for (l=0; l<max2; l++,ct2++) { 2264 #if defined(PETSC_USE_CTABLE) 2265 ierr = PetscTableFind(cmap_i,rbuf3_i[ct2]+1,&tcol);CHKERRQ(ierr); 2266 if (tcol) { 2267 #else 2268 if ((tcol = cmap_i[rbuf3_i[ct2]])) { 2269 #endif 2270 *mat_j++ = tcol - 1; 2271 if (!ijonly) { 2272 ierr = PetscMemcpy(mat_a,rbuf4_i+ct2*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 2273 mat_a += bs2; 2274 } 2275 ilen++; 2276 } 2277 } 2278 } else { /* allcolumns */ 2279 for (l=0; l<max2; l++,ct2++) { 2280 *mat_j++ = rbuf3_i[ct2]; 2281 if (!ijonly) { 2282 ierr = PetscMemcpy(mat_a,rbuf4_i+ct2*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 2283 mat_a += bs2; 2284 } 2285 ilen++; 2286 } 2287 } 2288 imat_ilen[row] = ilen; 2289 } 2290 } 2291 } 2292 } 2293 /* sort the rows */ 2294 { 2295 PetscInt ilen_row,*imat_ilen,*imat_j,*imat_i; 2296 MatScalar *imat_a = NULL; 2297 MatScalar *work; 2298 2299 ierr = PetscMalloc1(bs2,&work);CHKERRQ(ierr); 2300 for (i=0; i<ismax; i++) { 2301 mat = (Mat_SeqBAIJ*)submats[i]->data; 2302 imat_ilen = mat->ilen; 2303 imat_j = mat->j; 2304 imat_i = mat->i; 2305 if (!ijonly) imat_a = mat->a; 2306 if (allcolumns[i]) continue; 2307 jmax = nrow[i]; 2308 for (j=0; j<jmax; j++) { 2309 mat_i = imat_i[j]; 2310 if (!ijonly) mat_a = imat_a + mat_i*bs2; 2311 mat_j = imat_j + mat_i; 2312 ilen_row = imat_ilen[j]; 2313 if (!ijonly) {ierr = PetscSortIntWithDataArray(ilen_row,mat_j,mat_a,bs2*sizeof(MatScalar),work);CHKERRQ(ierr);} 2314 else {ierr = PetscSortInt(ilen_row,mat_j);CHKERRQ(ierr);} 2315 } 2316 } 2317 ierr = PetscFree(work);CHKERRQ(ierr); 2318 } 2319 if (!ijonly) { 2320 ierr = PetscFree(r_status4);CHKERRQ(ierr); 2321 ierr = PetscFree(r_waits4);CHKERRQ(ierr); 2322 if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits4,s_status4);CHKERRQ(ierr);} 2323 ierr = PetscFree(s_waits4);CHKERRQ(ierr); 2324 ierr = PetscFree(s_status4);CHKERRQ(ierr); 2325 } 2326 2327 /* Restore the indices */ 2328 for (i=0; i<ismax; i++) { 2329 if (!allrows[i]) { 2330 ierr = ISRestoreIndices(isrow[i],irow+i);CHKERRQ(ierr); 2331 } 2332 if (!allcolumns[i]) { 2333 ierr = ISRestoreIndices(iscol[i],icol+i);CHKERRQ(ierr); 2334 } 2335 } 2336 2337 /* Destroy allocated memory */ 2338 #if defined(PETSC_USE_CTABLE) 2339 ierr = PetscFree4(irow,icol,nrow,ncol);CHKERRQ(ierr); 2340 #else 2341 ierr = PetscFree5(irow,icol,nrow,ncol,rtable);CHKERRQ(ierr); 2342 #endif 2343 ierr = PetscFree4(w1,w2,w3,w4);CHKERRQ(ierr); 2344 ierr = PetscFree(pa);CHKERRQ(ierr); 2345 2346 ierr = PetscFree4(sbuf1,ptr,tmp,ctr);CHKERRQ(ierr); 2347 ierr = PetscFree(sbuf1);CHKERRQ(ierr); 2348 ierr = PetscFree(rbuf2);CHKERRQ(ierr); 2349 for (i=0; i<nrqr; ++i) { 2350 ierr = PetscFree(sbuf2[i]);CHKERRQ(ierr); 2351 } 2352 for (i=0; i<nrqs; ++i) { 2353 ierr = PetscFree(rbuf3[i]);CHKERRQ(ierr); 2354 } 2355 ierr = PetscFree3(sbuf2,req_size,req_source);CHKERRQ(ierr); 2356 ierr = PetscFree(rbuf3);CHKERRQ(ierr); 2357 ierr = PetscFree(sbuf_aj[0]);CHKERRQ(ierr); 2358 ierr = PetscFree(sbuf_aj);CHKERRQ(ierr); 2359 if (!ijonly) { 2360 for (i=0; i<nrqs; ++i) {ierr = PetscFree(rbuf4[i]);CHKERRQ(ierr);} 2361 ierr = PetscFree(rbuf4);CHKERRQ(ierr); 2362 ierr = PetscFree(sbuf_aa[0]);CHKERRQ(ierr); 2363 ierr = PetscFree(sbuf_aa);CHKERRQ(ierr); 2364 } 2365 2366 #if defined(PETSC_USE_CTABLE) 2367 for (i=0; i<ismax; i++) { 2368 ierr = PetscTableDestroy((PetscTable*)&rmap[i]);CHKERRQ(ierr); 2369 } 2370 #endif 2371 ierr = PetscFree(rmap);CHKERRQ(ierr); 2372 2373 for (i=0; i<ismax; i++) { 2374 if (!allcolumns[i]) { 2375 #if defined(PETSC_USE_CTABLE) 2376 ierr = PetscTableDestroy((PetscTable*)&cmap[i]);CHKERRQ(ierr); 2377 #else 2378 ierr = PetscFree(cmap[i]);CHKERRQ(ierr); 2379 #endif 2380 } 2381 } 2382 ierr = PetscFree(cmap);CHKERRQ(ierr); 2383 ierr = PetscFree(lens);CHKERRQ(ierr); 2384 2385 for (i=0; i<ismax; i++) { 2386 ierr = MatAssemblyBegin(submats[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2387 ierr = MatAssemblyEnd(submats[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2388 } 2389 2390 c->ijonly = PETSC_FALSE; /* set back to the default */ 2391 PetscFunctionReturn(0); 2392 } 2393 2394