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,N=C->cmap->N,bs=C->rmap->bs; 524 Mat_SeqBAIJ *subc; 525 Mat_SubMat *smat; 526 527 PetscFunctionBegin; 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 /* Determine the number of stages through which submatrices are done */ 535 nmax = 20*1000000 / (c->Nbs * sizeof(PetscInt)); 536 if (!nmax) nmax = 1; 537 538 if (scall == MAT_INITIAL_MATRIX) { 539 nstages_local = ismax/nmax + ((ismax % nmax) ? 1 : 0); 540 541 /* Make sure every processor loops through the nstages */ 542 ierr = MPIU_Allreduce(&nstages_local,&nstages,1,MPIU_INT,MPI_MAX,PetscObjectComm((PetscObject)C));CHKERRQ(ierr); 543 544 ierr = PetscCalloc1(ismax+1,submat);CHKERRQ(ierr); 545 } else { 546 subc = (Mat_SeqBAIJ*)((*submat)[0]->data); 547 smat = subc->submatis1; 548 if (!smat) { 549 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"MatGetSubMatrices(...,MAT_REUSE_MATRIX,...) requires submat"); 550 } 551 nstages = smat->nstages; 552 } 553 554 for (i=0,pos=0; i<nstages; i++) { 555 if (pos+nmax <= ismax) max_no = nmax; 556 else if (pos == ismax) max_no = 0; 557 else max_no = ismax-pos; 558 559 ierr = MatGetSubMatrices_MPIBAIJ_local(C,max_no,isrow_block+pos,iscol_block+pos,scall,*submat+pos);CHKERRQ(ierr); 560 pos += max_no; 561 } 562 563 if (scall == MAT_INITIAL_MATRIX && ismax) { 564 /* save nstages for reuse */ 565 subc = (Mat_SeqBAIJ*)((*submat)[0]->data); 566 smat = subc->submatis1; 567 if (!smat) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"smat does not exit"); 568 smat->nstages = nstages; 569 } 570 571 for (i=0; i<ismax; i++) { 572 ierr = ISDestroy(&isrow_block[i]);CHKERRQ(ierr); 573 ierr = ISDestroy(&iscol_block[i]);CHKERRQ(ierr); 574 } 575 ierr = PetscFree2(isrow_block,iscol_block);CHKERRQ(ierr); 576 PetscFunctionReturn(0); 577 } 578 579 #if defined(PETSC_USE_CTABLE) 580 PetscErrorCode PetscGetProc(const PetscInt row, const PetscMPIInt size, const PetscInt proc_gnode[], PetscMPIInt *rank) 581 { 582 PetscInt nGlobalNd = proc_gnode[size]; 583 PetscMPIInt fproc; 584 PetscErrorCode ierr; 585 586 PetscFunctionBegin; 587 ierr = PetscMPIIntCast((PetscInt)(((float)row * (float)size / (float)nGlobalNd + 0.5)),&fproc);CHKERRQ(ierr); 588 if (fproc > size) fproc = size; 589 while (row < proc_gnode[fproc] || row >= proc_gnode[fproc+1]) { 590 if (row < proc_gnode[fproc]) fproc--; 591 else fproc++; 592 } 593 *rank = fproc; 594 PetscFunctionReturn(0); 595 } 596 #endif 597 598 /* -------------------------------------------------------------------------*/ 599 /* This code is used for BAIJ and SBAIJ matrices (unfortunate dependency) */ 600 PetscErrorCode MatGetSubMatrices_MPIBAIJ_local(Mat C,PetscInt ismax,const IS isrow[],const IS iscol[],MatReuse scall,Mat *submats) 601 { 602 Mat_MPIBAIJ *c = (Mat_MPIBAIJ*)C->data; 603 Mat A = c->A; 604 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)c->B->data,*subc; 605 const PetscInt **icol,**irow; 606 PetscInt *nrow,*ncol,start; 607 PetscErrorCode ierr; 608 PetscMPIInt rank,size,tag0,tag2,tag3,tag4,*w1,*w2,*w3,*w4,nrqr; 609 PetscInt **sbuf1,**sbuf2,*sbuf2_i,i,j,k,l,ct1,ct2,**rbuf1,row,proc=-1; 610 PetscInt nrqs=0,msz,**ptr=NULL,*req_size=NULL,*ctr=NULL,*pa,*tmp=NULL,tcol; 611 PetscInt **rbuf3=NULL,*req_source1=NULL,*req_source2,**sbuf_aj,**rbuf2=NULL,max1,max2; 612 PetscInt **lens,is_no,ncols,*cols,mat_i,*mat_j,tmp2,jmax; 613 #if defined(PETSC_USE_CTABLE) 614 PetscTable *cmap,cmap_i=NULL,*rmap,rmap_i; 615 #else 616 PetscInt **cmap,*cmap_i=NULL,**rmap,*rmap_i; 617 #endif 618 const PetscInt *irow_i,*icol_i; 619 PetscInt ctr_j,*sbuf1_j,*sbuf_aj_i,*rbuf1_i,kmax,*lens_i; 620 MPI_Request *s_waits1,*r_waits1,*s_waits2,*r_waits2,*r_waits3; 621 MPI_Request *r_waits4,*s_waits3,*s_waits4; 622 MPI_Status *r_status1,*r_status2,*s_status1,*s_status3,*s_status2; 623 MPI_Status *r_status3,*r_status4,*s_status4; 624 MPI_Comm comm; 625 PetscScalar **rbuf4,*rbuf4_i=NULL,**sbuf_aa,*vals,*mat_a=NULL,*imat_a,*sbuf_aa_i; 626 PetscMPIInt *onodes1,*olengths1,end; 627 PetscInt **row2proc,*row2proc_i,*imat_ilen,*imat_j,*imat_i; 628 Mat_SubMat **smats,*smat_i; 629 PetscBool *issorted,colflag,iscsorted=PETSC_TRUE; 630 PetscInt *sbuf1_i,*rbuf2_i,*rbuf3_i,ilen; 631 PetscInt bs=C->rmap->bs,bs2=c->bs2,rstart = c->rstartbs; 632 PetscBool ijonly=c->ijonly; /* private flag indicates only matrix data structures are requested */ 633 PetscInt nzA,nzB,*a_i=a->i,*b_i=b->i,*a_j = a->j,*b_j = b->j,ctmp,imark,*cworkA,*cworkB; 634 PetscScalar *vworkA=NULL,*vworkB=NULL,*a_a = a->a,*b_a = b->a; 635 PetscInt cstart = c->cstartbs,*bmap = c->garray; 636 PetscBool *allrows,*allcolumns; 637 638 PetscFunctionBegin; 639 ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr); 640 size = c->size; 641 rank = c->rank; 642 643 ierr = PetscMalloc6(ismax,&smats,ismax,&row2proc,ismax,&cmap,ismax,&rmap,ismax,&allcolumns,ismax,&allrows);CHKERRQ(ierr); 644 ierr = PetscMalloc5(ismax,&irow,ismax,&icol,ismax,&nrow,ismax,&ncol,ismax,&issorted);CHKERRQ(ierr); 645 646 for (i=0; i<ismax; i++) { 647 ierr = ISSorted(iscol[i],&issorted[i]);CHKERRQ(ierr); 648 if (!issorted[i]) iscsorted = issorted[i]; /* columns are not sorted! */ 649 ierr = ISSorted(isrow[i],&issorted[i]);CHKERRQ(ierr); 650 651 /* Check for special case: allcolumns */ 652 ierr = ISIdentity(iscol[i],&colflag);CHKERRQ(ierr); 653 ierr = ISGetLocalSize(iscol[i],&ncol[i]);CHKERRQ(ierr); 654 655 if (colflag && ncol[i] == c->Nbs) { 656 allcolumns[i] = PETSC_TRUE; 657 icol[i] = NULL; 658 } else { 659 allcolumns[i] = PETSC_FALSE; 660 ierr = ISGetIndices(iscol[i],&icol[i]);CHKERRQ(ierr); 661 } 662 663 /* Check for special case: allrows */ 664 ierr = ISIdentity(isrow[i],&colflag);CHKERRQ(ierr); 665 ierr = ISGetLocalSize(isrow[i],&nrow[i]);CHKERRQ(ierr); 666 if (colflag && nrow[i] == c->Mbs) { 667 allrows[i] = PETSC_TRUE; 668 irow[i] = NULL; 669 } else { 670 allrows[i] = PETSC_FALSE; 671 ierr = ISGetIndices(isrow[i],&irow[i]);CHKERRQ(ierr); 672 } 673 } 674 675 if (scall == MAT_REUSE_MATRIX) { 676 /* Assumes new rows are same length as the old rows */ 677 for (i=0; i<ismax; i++) { 678 subc = (Mat_SeqBAIJ*)(submats[i]->data); 679 if (subc->mbs != nrow[i] || subc->nbs != ncol[i]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong size"); 680 681 /* Initial matrix as if empty */ 682 ierr = PetscMemzero(subc->ilen,subc->mbs*sizeof(PetscInt));CHKERRQ(ierr); 683 684 /* Initial matrix as if empty */ 685 submats[i]->factortype = C->factortype; 686 687 smat_i = subc->submatis1; 688 smats[i] = smat_i; 689 690 nrqs = smat_i->nrqs; 691 nrqr = smat_i->nrqr; 692 rbuf1 = smat_i->rbuf1; 693 rbuf2 = smat_i->rbuf2; 694 rbuf3 = smat_i->rbuf3; 695 req_source2 = smat_i->req_source2; 696 697 sbuf1 = smat_i->sbuf1; 698 sbuf2 = smat_i->sbuf2; 699 ptr = smat_i->ptr; 700 tmp = smat_i->tmp; 701 ctr = smat_i->ctr; 702 703 pa = smat_i->pa; 704 req_size = smat_i->req_size; 705 req_source1 = smat_i->req_source1; 706 707 allcolumns[i] = smat_i->allcolumns; 708 allrows[i] = smat_i->allrows; 709 row2proc[i] = smat_i->row2proc; 710 rmap[i] = smat_i->rmap; 711 cmap[i] = smat_i->cmap; 712 } 713 } else { /* scall == MAT_INITIAL_MATRIX */ 714 /* Get some new tags to keep the communication clean */ 715 ierr = PetscObjectGetNewTag((PetscObject)C,&tag2);CHKERRQ(ierr); 716 ierr = PetscObjectGetNewTag((PetscObject)C,&tag3);CHKERRQ(ierr); 717 718 /* evaluate communication - mesg to who, length of mesg, and buffer space 719 required. Based on this, buffers are allocated, and data copied into them*/ 720 ierr = PetscCalloc4(size,&w1,size,&w2,size,&w3,size,&w4);CHKERRQ(ierr); /* mesg size, initialize work vectors */ 721 722 for (i=0; i<ismax; i++) { 723 jmax = nrow[i]; 724 irow_i = irow[i]; 725 726 ierr = PetscMalloc1(jmax,&row2proc_i);CHKERRQ(ierr); 727 row2proc[i] = row2proc_i; 728 729 if (issorted[i]) proc = 0; 730 for (j=0; j<jmax; j++) { 731 if (!issorted[i]) proc = 0; 732 if (allrows[i]) row = j; 733 else row = irow_i[j]; 734 735 while (row >= c->rangebs[proc+1]) proc++; 736 w4[proc]++; 737 row2proc_i[j] = proc; /* map row index to proc */ 738 } 739 for (j=0; j<size; j++) { 740 if (w4[j]) { w1[j] += w4[j]; w3[j]++; w4[j] = 0;} 741 } 742 } 743 744 nrqs = 0; /* no of outgoing messages */ 745 msz = 0; /* total mesg length (for all procs) */ 746 w1[rank] = 0; /* no mesg sent to self */ 747 w3[rank] = 0; 748 for (i=0; i<size; i++) { 749 if (w1[i]) { w2[i] = 1; nrqs++;} /* there exists a message to proc i */ 750 } 751 ierr = PetscMalloc1(nrqs+1,&pa);CHKERRQ(ierr); /*(proc -array)*/ 752 for (i=0,j=0; i<size; i++) { 753 if (w1[i]) { pa[j] = i; j++; } 754 } 755 756 /* Each message would have a header = 1 + 2*(no of IS) + data */ 757 for (i=0; i<nrqs; i++) { 758 j = pa[i]; 759 w1[j] += w2[j] + 2* w3[j]; 760 msz += w1[j]; 761 } 762 ierr = PetscInfo2(0,"Number of outgoing messages %D Total message length %D\n",nrqs,msz);CHKERRQ(ierr); 763 764 /* Determine the number of messages to expect, their lengths, from from-ids */ 765 ierr = PetscGatherNumberOfMessages(comm,w2,w1,&nrqr);CHKERRQ(ierr); 766 ierr = PetscGatherMessageLengths(comm,nrqs,nrqr,w1,&onodes1,&olengths1);CHKERRQ(ierr); 767 768 /* Now post the Irecvs corresponding to these messages */ 769 tag0 = ((PetscObject)C)->tag; 770 ierr = PetscPostIrecvInt(comm,tag0,nrqr,onodes1,olengths1,&rbuf1,&r_waits1);CHKERRQ(ierr); 771 772 ierr = PetscFree(onodes1);CHKERRQ(ierr); 773 ierr = PetscFree(olengths1);CHKERRQ(ierr); 774 775 /* Allocate Memory for outgoing messages */ 776 ierr = PetscMalloc4(size,&sbuf1,size,&ptr,2*msz,&tmp,size,&ctr);CHKERRQ(ierr); 777 ierr = PetscMemzero(sbuf1,size*sizeof(PetscInt*));CHKERRQ(ierr); 778 ierr = PetscMemzero(ptr,size*sizeof(PetscInt*));CHKERRQ(ierr); 779 780 { 781 PetscInt *iptr = tmp; 782 k = 0; 783 for (i=0; i<nrqs; i++) { 784 j = pa[i]; 785 iptr += k; 786 sbuf1[j] = iptr; 787 k = w1[j]; 788 } 789 } 790 791 /* Form the outgoing messages. Initialize the header space */ 792 for (i=0; i<nrqs; i++) { 793 j = pa[i]; 794 sbuf1[j][0] = 0; 795 ierr = PetscMemzero(sbuf1[j]+1,2*w3[j]*sizeof(PetscInt));CHKERRQ(ierr); 796 ptr[j] = sbuf1[j] + 2*w3[j] + 1; 797 } 798 799 /* Parse the isrow and copy data into outbuf */ 800 for (i=0; i<ismax; i++) { 801 row2proc_i = row2proc[i]; 802 ierr = PetscMemzero(ctr,size*sizeof(PetscInt));CHKERRQ(ierr); 803 irow_i = irow[i]; 804 jmax = nrow[i]; 805 for (j=0; j<jmax; j++) { /* parse the indices of each IS */ 806 proc = row2proc_i[j]; 807 if (allrows[i]) row = j; 808 else row = irow_i[j]; 809 810 if (proc != rank) { /* copy to the outgoing buf*/ 811 ctr[proc]++; 812 *ptr[proc] = row; 813 ptr[proc]++; 814 } 815 } 816 /* Update the headers for the current IS */ 817 for (j=0; j<size; j++) { /* Can Optimise this loop too */ 818 if ((ctr_j = ctr[j])) { 819 sbuf1_j = sbuf1[j]; 820 k = ++sbuf1_j[0]; 821 sbuf1_j[2*k] = ctr_j; 822 sbuf1_j[2*k-1] = i; 823 } 824 } 825 } 826 827 /* Now post the sends */ 828 ierr = PetscMalloc1(nrqs+1,&s_waits1);CHKERRQ(ierr); 829 for (i=0; i<nrqs; ++i) { 830 j = pa[i]; 831 ierr = MPI_Isend(sbuf1[j],w1[j],MPIU_INT,j,tag0,comm,s_waits1+i);CHKERRQ(ierr); 832 } 833 834 /* Post Receives to capture the buffer size */ 835 ierr = PetscMalloc1(nrqs+1,&r_waits2);CHKERRQ(ierr); 836 ierr = PetscMalloc3(nrqs+1,&req_source2,nrqs+1,&rbuf2,nrqs+1,&rbuf3);CHKERRQ(ierr); 837 rbuf2[0] = tmp + msz; 838 for (i=1; i<nrqs; ++i) { 839 rbuf2[i] = rbuf2[i-1]+w1[pa[i-1]]; 840 } 841 for (i=0; i<nrqs; ++i) { 842 j = pa[i]; 843 ierr = MPI_Irecv(rbuf2[i],w1[j],MPIU_INT,j,tag2,comm,r_waits2+i);CHKERRQ(ierr); 844 } 845 846 /* Send to other procs the buf size they should allocate */ 847 /* Receive messages*/ 848 ierr = PetscMalloc1(nrqr+1,&s_waits2);CHKERRQ(ierr); 849 ierr = PetscMalloc1(nrqr+1,&r_status1);CHKERRQ(ierr); 850 ierr = PetscMalloc3(nrqr,&sbuf2,nrqr,&req_size,nrqr,&req_source1);CHKERRQ(ierr); 851 852 ierr = MPI_Waitall(nrqr,r_waits1,r_status1);CHKERRQ(ierr); 853 for (i=0; i<nrqr; ++i) { 854 req_size[i] = 0; 855 rbuf1_i = rbuf1[i]; 856 start = 2*rbuf1_i[0] + 1; 857 ierr = MPI_Get_count(r_status1+i,MPIU_INT,&end);CHKERRQ(ierr); 858 ierr = PetscMalloc1(end+1,&sbuf2[i]);CHKERRQ(ierr); 859 sbuf2_i = sbuf2[i]; 860 for (j=start; j<end; j++) { 861 row = rbuf1_i[j] - rstart; 862 ncols = a_i[row+1] - a_i[row] + b_i[row+1] - b_i[row]; 863 sbuf2_i[j] = ncols; 864 req_size[i] += ncols; 865 } 866 req_source1[i] = r_status1[i].MPI_SOURCE; 867 /* form the header */ 868 sbuf2_i[0] = req_size[i]; 869 for (j=1; j<start; j++) sbuf2_i[j] = rbuf1_i[j]; 870 871 ierr = MPI_Isend(sbuf2_i,end,MPIU_INT,req_source1[i],tag2,comm,s_waits2+i);CHKERRQ(ierr); 872 } 873 874 ierr = PetscFree(r_status1);CHKERRQ(ierr); 875 ierr = PetscFree(r_waits1);CHKERRQ(ierr); 876 ierr = PetscFree4(w1,w2,w3,w4);CHKERRQ(ierr); 877 878 /* Receive messages*/ 879 ierr = PetscMalloc1(nrqs+1,&r_waits3);CHKERRQ(ierr); 880 ierr = PetscMalloc1(nrqs+1,&r_status2);CHKERRQ(ierr); 881 882 ierr = MPI_Waitall(nrqs,r_waits2,r_status2);CHKERRQ(ierr); 883 for (i=0; i<nrqs; ++i) { 884 ierr = PetscMalloc1(rbuf2[i][0]+1,&rbuf3[i]);CHKERRQ(ierr); 885 req_source2[i] = r_status2[i].MPI_SOURCE; 886 ierr = MPI_Irecv(rbuf3[i],rbuf2[i][0],MPIU_INT,req_source2[i],tag3,comm,r_waits3+i);CHKERRQ(ierr); 887 } 888 ierr = PetscFree(r_status2);CHKERRQ(ierr); 889 ierr = PetscFree(r_waits2);CHKERRQ(ierr); 890 891 /* Wait on sends1 and sends2 */ 892 ierr = PetscMalloc1(nrqs+1,&s_status1);CHKERRQ(ierr); 893 ierr = PetscMalloc1(nrqr+1,&s_status2);CHKERRQ(ierr); 894 895 if (nrqs) {ierr = MPI_Waitall(nrqs,s_waits1,s_status1);CHKERRQ(ierr);} 896 if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits2,s_status2);CHKERRQ(ierr);} 897 ierr = PetscFree(s_status1);CHKERRQ(ierr); 898 ierr = PetscFree(s_status2);CHKERRQ(ierr); 899 ierr = PetscFree(s_waits1);CHKERRQ(ierr); 900 ierr = PetscFree(s_waits2);CHKERRQ(ierr); 901 902 /* Now allocate sending buffers for a->j, and send them off */ 903 ierr = PetscMalloc1(nrqr+1,&sbuf_aj);CHKERRQ(ierr); 904 for (i=0,j=0; i<nrqr; i++) j += req_size[i]; 905 ierr = PetscMalloc1(j+1,&sbuf_aj[0]);CHKERRQ(ierr); 906 for (i=1; i<nrqr; i++) sbuf_aj[i] = sbuf_aj[i-1] + req_size[i-1]; 907 908 ierr = PetscMalloc1(nrqr+1,&s_waits3);CHKERRQ(ierr); 909 { 910 911 for (i=0; i<nrqr; i++) { 912 rbuf1_i = rbuf1[i]; 913 sbuf_aj_i = sbuf_aj[i]; 914 ct1 = 2*rbuf1_i[0] + 1; 915 ct2 = 0; 916 for (j=1,max1=rbuf1_i[0]; j<=max1; j++) { 917 kmax = rbuf1[i][2*j]; 918 for (k=0; k<kmax; k++,ct1++) { 919 row = rbuf1_i[ct1] - rstart; 920 nzA = a_i[row+1] - a_i[row]; nzB = b_i[row+1] - b_i[row]; 921 ncols = nzA + nzB; 922 cworkA = a_j + a_i[row]; cworkB = b_j + b_i[row]; 923 924 /* load the column indices for this row into cols */ 925 cols = sbuf_aj_i + ct2; 926 for (l=0; l<nzB; l++) { 927 if ((ctmp = bmap[cworkB[l]]) < cstart) cols[l] = ctmp; 928 else break; 929 } 930 imark = l; 931 for (l=0; l<nzA; l++) {cols[imark+l] = cstart + cworkA[l];} 932 for (l=imark; l<nzB; l++) cols[nzA+l] = bmap[cworkB[l]]; 933 ct2 += ncols; 934 } 935 } 936 ierr = MPI_Isend(sbuf_aj_i,req_size[i],MPIU_INT,req_source1[i],tag3,comm,s_waits3+i);CHKERRQ(ierr); 937 } 938 } 939 ierr = PetscMalloc2(nrqs+1,&r_status3,nrqr+1,&s_status3);CHKERRQ(ierr); 940 941 /* create col map: global col of C -> local col of submatrices */ 942 #if defined(PETSC_USE_CTABLE) 943 for (i=0; i<ismax; i++) { 944 if (!allcolumns[i]) { 945 ierr = PetscTableCreate(ncol[i]+1,c->Nbs+1,&cmap[i]);CHKERRQ(ierr); 946 947 jmax = ncol[i]; 948 icol_i = icol[i]; 949 cmap_i = cmap[i]; 950 for (j=0; j<jmax; j++) { 951 ierr = PetscTableAdd(cmap[i],icol_i[j]+1,j+1,INSERT_VALUES);CHKERRQ(ierr); 952 } 953 } else cmap[i] = NULL; 954 } 955 #else 956 for (i=0; i<ismax; i++) { 957 if (!allcolumns[i]) { 958 ierr = PetscCalloc1(c->Nbs,&cmap[i]);CHKERRQ(ierr); 959 jmax = ncol[i]; 960 icol_i = icol[i]; 961 cmap_i = cmap[i]; 962 for (j=0; j<jmax; j++) cmap_i[icol_i[j]] = j+1; 963 } else cmap[i] = NULL; 964 } 965 #endif 966 967 /* Create lens which is required for MatCreate... */ 968 for (i=0,j=0; i<ismax; i++) j += nrow[i]; 969 ierr = PetscMalloc1(ismax,&lens);CHKERRQ(ierr); 970 971 if (ismax) { 972 ierr = PetscCalloc1(j,&lens[0]);CHKERRQ(ierr); 973 } 974 for (i=1; i<ismax; i++) lens[i] = lens[i-1] + nrow[i-1]; 975 976 /* Update lens from local data */ 977 for (i=0; i<ismax; i++) { 978 row2proc_i = row2proc[i]; 979 jmax = nrow[i]; 980 if (!allcolumns[i]) cmap_i = cmap[i]; 981 irow_i = irow[i]; 982 lens_i = lens[i]; 983 for (j=0; j<jmax; j++) { 984 if (allrows[i]) row = j; 985 else row = irow_i[j]; /* global blocked row of C */ 986 987 proc = row2proc_i[j]; 988 if (proc == rank) { 989 /* Get indices from matA and then from matB */ 990 #if defined(PETSC_USE_CTABLE) 991 PetscInt tt; 992 #endif 993 row = row - rstart; 994 nzA = a_i[row+1] - a_i[row]; 995 nzB = b_i[row+1] - b_i[row]; 996 cworkA = a_j + a_i[row]; 997 cworkB = b_j + b_i[row]; 998 999 if (!allcolumns[i]) { 1000 #if defined(PETSC_USE_CTABLE) 1001 for (k=0; k<nzA; k++) { 1002 ierr = PetscTableFind(cmap_i,cstart+cworkA[k]+1,&tt);CHKERRQ(ierr); 1003 if (tt) lens_i[j]++; 1004 } 1005 for (k=0; k<nzB; k++) { 1006 ierr = PetscTableFind(cmap_i,bmap[cworkB[k]]+1,&tt);CHKERRQ(ierr); 1007 if (tt) lens_i[j]++; 1008 } 1009 1010 #else 1011 for (k=0; k<nzA; k++) { 1012 if (cmap_i[cstart + cworkA[k]]) lens_i[j]++; 1013 } 1014 for (k=0; k<nzB; k++) { 1015 if (cmap_i[bmap[cworkB[k]]]) lens_i[j]++; 1016 } 1017 #endif 1018 } else { /* allcolumns */ 1019 lens_i[j] = nzA + nzB; 1020 } 1021 } 1022 } 1023 } 1024 1025 /* Create row map: global row of C -> local row of submatrices */ 1026 for (i=0; i<ismax; i++) { 1027 if (!allrows[i]) { 1028 #if defined(PETSC_USE_CTABLE) 1029 ierr = PetscTableCreate(nrow[i]+1,c->Mbs+1,&rmap[i]);CHKERRQ(ierr); 1030 irow_i = irow[i]; 1031 jmax = nrow[i]; 1032 for (j=0; j<jmax; j++) { 1033 if (allrows[i]) { 1034 ierr = PetscTableAdd(rmap[i],j+1,j+1,INSERT_VALUES);CHKERRQ(ierr); 1035 } else { 1036 ierr = PetscTableAdd(rmap[i],irow_i[j]+1,j+1,INSERT_VALUES);CHKERRQ(ierr); 1037 } 1038 } 1039 #else 1040 ierr = PetscCalloc1(c->Mbs,&rmap[i]);CHKERRQ(ierr); 1041 rmap_i = rmap[i]; 1042 irow_i = irow[i]; 1043 jmax = nrow[i]; 1044 for (j=0; j<jmax; j++) { 1045 if (allrows[i]) rmap_i[j] = j; 1046 else rmap_i[irow_i[j]] = j; 1047 } 1048 #endif 1049 } else rmap[i] = NULL; 1050 } 1051 1052 /* Update lens from offproc data */ 1053 { 1054 PetscInt *rbuf2_i,*rbuf3_i,*sbuf1_i; 1055 1056 ierr = MPI_Waitall(nrqs,r_waits3,r_status3);CHKERRQ(ierr); 1057 for (tmp2=0; tmp2<nrqs; tmp2++) { 1058 sbuf1_i = sbuf1[pa[tmp2]]; 1059 jmax = sbuf1_i[0]; 1060 ct1 = 2*jmax+1; 1061 ct2 = 0; 1062 rbuf2_i = rbuf2[tmp2]; 1063 rbuf3_i = rbuf3[tmp2]; 1064 for (j=1; j<=jmax; j++) { 1065 is_no = sbuf1_i[2*j-1]; 1066 max1 = sbuf1_i[2*j]; 1067 lens_i = lens[is_no]; 1068 if (!allcolumns[is_no]) cmap_i = cmap[is_no]; 1069 rmap_i = rmap[is_no]; 1070 for (k=0; k<max1; k++,ct1++) { 1071 if (allrows[is_no]) { 1072 row = sbuf1_i[ct1]; 1073 } else { 1074 #if defined(PETSC_USE_CTABLE) 1075 ierr = PetscTableFind(rmap_i,sbuf1_i[ct1]+1,&row);CHKERRQ(ierr); 1076 row--; 1077 if (row < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"row not found in table"); 1078 #else 1079 row = rmap_i[sbuf1_i[ct1]]; /* the val in the new matrix to be */ 1080 #endif 1081 } 1082 max2 = rbuf2_i[ct1]; 1083 for (l=0; l<max2; l++,ct2++) { 1084 if (!allcolumns[is_no]) { 1085 #if defined(PETSC_USE_CTABLE) 1086 ierr = PetscTableFind(cmap_i,rbuf3_i[ct2]+1,&tcol);CHKERRQ(ierr); 1087 #else 1088 tcol = cmap_i[rbuf3_i[ct2]]; 1089 #endif 1090 if (tcol) lens_i[row]++; 1091 } else { /* allcolumns */ 1092 lens_i[row]++; /* lens_i[row] += max2 ? */ 1093 } 1094 } 1095 } 1096 } 1097 } 1098 } 1099 ierr = PetscFree(r_waits3);CHKERRQ(ierr); 1100 if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits3,s_status3);CHKERRQ(ierr);} 1101 ierr = PetscFree2(r_status3,s_status3);CHKERRQ(ierr); 1102 ierr = PetscFree(s_waits3);CHKERRQ(ierr); 1103 1104 /* Create the submatrices */ 1105 for (i=0; i<ismax; i++) { 1106 PetscInt bs_tmp; 1107 if (ijonly) bs_tmp = 1; 1108 else bs_tmp = bs; 1109 1110 ierr = MatCreate(PETSC_COMM_SELF,submats+i);CHKERRQ(ierr); 1111 ierr = MatSetSizes(submats[i],nrow[i]*bs_tmp,ncol[i]*bs_tmp,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 1112 1113 ierr = MatSetType(submats[i],((PetscObject)A)->type_name);CHKERRQ(ierr); 1114 ierr = MatSeqBAIJSetPreallocation(submats[i],bs_tmp,0,lens[i]);CHKERRQ(ierr); 1115 ierr = MatSeqSBAIJSetPreallocation(submats[i],bs_tmp,0,lens[i]);CHKERRQ(ierr); /* this subroutine is used by SBAIJ routines */ 1116 1117 /* create struct Mat_SubMat and attached it to submat */ 1118 ierr = PetscNew(&smat_i);CHKERRQ(ierr); 1119 subc = (Mat_SeqBAIJ*)submats[i]->data; 1120 subc->submatis1 = smat_i; 1121 smats[i] = smat_i; 1122 1123 smat_i->destroy = submats[i]->ops->destroy; 1124 submats[i]->ops->destroy = MatDestroy_SeqBAIJ_Submatrices; 1125 submats[i]->factortype = C->factortype; 1126 1127 smat_i->id = i; 1128 smat_i->nrqs = nrqs; 1129 smat_i->nrqr = nrqr; 1130 smat_i->rbuf1 = rbuf1; 1131 smat_i->rbuf2 = rbuf2; 1132 smat_i->rbuf3 = rbuf3; 1133 smat_i->sbuf2 = sbuf2; 1134 smat_i->req_source2 = req_source2; 1135 1136 smat_i->sbuf1 = sbuf1; 1137 smat_i->ptr = ptr; 1138 smat_i->tmp = tmp; 1139 smat_i->ctr = ctr; 1140 1141 smat_i->pa = pa; 1142 smat_i->req_size = req_size; 1143 smat_i->req_source1 = req_source1; 1144 1145 smat_i->allcolumns = allcolumns[i]; 1146 smat_i->allrows = allrows[i]; 1147 smat_i->singleis = PETSC_FALSE; 1148 smat_i->row2proc = row2proc[i]; 1149 smat_i->rmap = rmap[i]; 1150 smat_i->cmap = cmap[i]; 1151 } 1152 1153 if (ismax) {ierr = PetscFree(lens[0]);CHKERRQ(ierr);} 1154 ierr = PetscFree(lens);CHKERRQ(ierr); 1155 ierr = PetscFree(sbuf_aj[0]);CHKERRQ(ierr); 1156 ierr = PetscFree(sbuf_aj);CHKERRQ(ierr); 1157 1158 } /* endof scall == MAT_INITIAL_MATRIX */ 1159 1160 /* Post recv matrix values */ 1161 if (!ijonly) { 1162 ierr = PetscObjectGetNewTag((PetscObject)C,&tag4);CHKERRQ(ierr); 1163 ierr = PetscMalloc1(nrqs+1,&rbuf4);CHKERRQ(ierr); 1164 ierr = PetscMalloc1(nrqs+1,&r_waits4);CHKERRQ(ierr); 1165 ierr = PetscMalloc1(nrqs+1,&r_status4);CHKERRQ(ierr); 1166 ierr = PetscMalloc1(nrqr+1,&s_status4);CHKERRQ(ierr); 1167 for (i=0; i<nrqs; ++i) { 1168 ierr = PetscMalloc1(rbuf2[i][0]*bs2,&rbuf4[i]);CHKERRQ(ierr); 1169 ierr = MPI_Irecv(rbuf4[i],rbuf2[i][0]*bs2,MPIU_SCALAR,req_source2[i],tag4,comm,r_waits4+i);CHKERRQ(ierr); 1170 } 1171 1172 /* Allocate sending buffers for a->a, and send them off */ 1173 ierr = PetscMalloc1(nrqr+1,&sbuf_aa);CHKERRQ(ierr); 1174 for (i=0,j=0; i<nrqr; i++) j += req_size[i]; 1175 1176 ierr = PetscMalloc1((j+1)*bs2,&sbuf_aa[0]);CHKERRQ(ierr); 1177 for (i=1; i<nrqr; i++) sbuf_aa[i] = sbuf_aa[i-1] + req_size[i-1]*bs2; 1178 1179 ierr = PetscMalloc1(nrqr+1,&s_waits4);CHKERRQ(ierr); 1180 1181 for (i=0; i<nrqr; i++) { 1182 rbuf1_i = rbuf1[i]; 1183 sbuf_aa_i = sbuf_aa[i]; 1184 ct1 = 2*rbuf1_i[0]+1; 1185 ct2 = 0; 1186 for (j=1,max1=rbuf1_i[0]; j<=max1; j++) { 1187 kmax = rbuf1_i[2*j]; 1188 for (k=0; k<kmax; k++,ct1++) { 1189 row = rbuf1_i[ct1] - rstart; 1190 nzA = a_i[row+1] - a_i[row]; 1191 nzB = b_i[row+1] - b_i[row]; 1192 ncols = nzA + nzB; 1193 cworkB = b_j + b_i[row]; 1194 vworkA = a_a + a_i[row]*bs2; 1195 vworkB = b_a + b_i[row]*bs2; 1196 1197 /* load the column values for this row into vals*/ 1198 vals = sbuf_aa_i+ct2*bs2; 1199 for (l=0; l<nzB; l++) { 1200 if ((bmap[cworkB[l]]) < cstart) { 1201 ierr = PetscMemcpy(vals+l*bs2,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1202 } else break; 1203 } 1204 imark = l; 1205 for (l=0; l<nzA; l++) { 1206 ierr = PetscMemcpy(vals+(imark+l)*bs2,vworkA+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1207 } 1208 for (l=imark; l<nzB; l++) { 1209 ierr = PetscMemcpy(vals+(nzA+l)*bs2,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1210 } 1211 1212 ct2 += ncols; 1213 } 1214 } 1215 ierr = MPI_Isend(sbuf_aa_i,req_size[i]*bs2,MPIU_SCALAR,req_source1[i],tag4,comm,s_waits4+i);CHKERRQ(ierr); 1216 } 1217 } 1218 1219 if (!ismax) { 1220 ierr = PetscFree(rbuf1[0]);CHKERRQ(ierr); 1221 ierr = PetscFree(rbuf1);CHKERRQ(ierr); 1222 } 1223 1224 /* Assemble the matrices */ 1225 /* First assemble the local rows */ 1226 for (i=0; i<ismax; i++) { 1227 row2proc_i = row2proc[i]; 1228 subc = (Mat_SeqBAIJ*)submats[i]->data; 1229 imat_ilen = subc->ilen; 1230 imat_j = subc->j; 1231 imat_i = subc->i; 1232 imat_a = subc->a; 1233 1234 if (!allcolumns[i]) cmap_i = cmap[i]; 1235 rmap_i = rmap[i]; 1236 irow_i = irow[i]; 1237 jmax = nrow[i]; 1238 for (j=0; j<jmax; j++) { 1239 if (allrows[i]) row = j; 1240 else row = irow_i[j]; 1241 proc = row2proc_i[j]; 1242 1243 if (proc == rank) { 1244 1245 row = row - rstart; 1246 nzA = a_i[row+1] - a_i[row]; 1247 nzB = b_i[row+1] - b_i[row]; 1248 cworkA = a_j + a_i[row]; 1249 cworkB = b_j + b_i[row]; 1250 if (!ijonly) { 1251 vworkA = a_a + a_i[row]*bs2; 1252 vworkB = b_a + b_i[row]*bs2; 1253 } 1254 1255 if (allrows[i]) { 1256 row = row+rstart; 1257 } else { 1258 #if defined(PETSC_USE_CTABLE) 1259 ierr = PetscTableFind(rmap_i,row+rstart+1,&row);CHKERRQ(ierr); 1260 row--; 1261 1262 if (row < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"row not found in table"); 1263 #else 1264 row = rmap_i[row + rstart]; 1265 #endif 1266 } 1267 mat_i = imat_i[row]; 1268 if (!ijonly) mat_a = imat_a + mat_i*bs2; 1269 mat_j = imat_j + mat_i; 1270 ilen = imat_ilen[row]; 1271 1272 /* load the column indices for this row into cols*/ 1273 if (!allcolumns[i]) { 1274 for (l=0; l<nzB; l++) { 1275 if ((ctmp = bmap[cworkB[l]]) < cstart) { 1276 #if defined(PETSC_USE_CTABLE) 1277 ierr = PetscTableFind(cmap_i,ctmp+1,&tcol);CHKERRQ(ierr); 1278 if (tcol) { 1279 #else 1280 if ((tcol = cmap_i[ctmp])) { 1281 #endif 1282 *mat_j++ = tcol - 1; 1283 ierr = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1284 mat_a += bs2; 1285 ilen++; 1286 } 1287 } else break; 1288 } 1289 imark = l; 1290 for (l=0; l<nzA; l++) { 1291 #if defined(PETSC_USE_CTABLE) 1292 ierr = PetscTableFind(cmap_i,cstart+cworkA[l]+1,&tcol);CHKERRQ(ierr); 1293 if (tcol) { 1294 #else 1295 if ((tcol = cmap_i[cstart + cworkA[l]])) { 1296 #endif 1297 *mat_j++ = tcol - 1; 1298 if (!ijonly) { 1299 ierr = PetscMemcpy(mat_a,vworkA+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1300 mat_a += bs2; 1301 } 1302 ilen++; 1303 } 1304 } 1305 for (l=imark; l<nzB; l++) { 1306 #if defined(PETSC_USE_CTABLE) 1307 ierr = PetscTableFind(cmap_i,bmap[cworkB[l]]+1,&tcol);CHKERRQ(ierr); 1308 if (tcol) { 1309 #else 1310 if ((tcol = cmap_i[bmap[cworkB[l]]])) { 1311 #endif 1312 *mat_j++ = tcol - 1; 1313 if (!ijonly) { 1314 ierr = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1315 mat_a += bs2; 1316 } 1317 ilen++; 1318 } 1319 } 1320 } else { /* allcolumns */ 1321 for (l=0; l<nzB; l++) { 1322 if ((ctmp = bmap[cworkB[l]]) < cstart) { 1323 *mat_j++ = ctmp; 1324 ierr = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1325 mat_a += bs2; 1326 ilen++; 1327 } else break; 1328 } 1329 imark = l; 1330 for (l=0; l<nzA; l++) { 1331 *mat_j++ = cstart+cworkA[l]; 1332 if (!ijonly) { 1333 ierr = PetscMemcpy(mat_a,vworkA+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1334 mat_a += bs2; 1335 } 1336 ilen++; 1337 } 1338 for (l=imark; l<nzB; l++) { 1339 *mat_j++ = bmap[cworkB[l]]; 1340 if (!ijonly) { 1341 ierr = PetscMemcpy(mat_a,vworkB+l*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1342 mat_a += bs2; 1343 } 1344 ilen++; 1345 } 1346 } 1347 imat_ilen[row] = ilen; 1348 } 1349 } 1350 } 1351 1352 /* Now assemble the off proc rows */ 1353 if (!ijonly) { 1354 ierr = MPI_Waitall(nrqs,r_waits4,r_status4);CHKERRQ(ierr); 1355 } 1356 for (tmp2=0; tmp2<nrqs; tmp2++) { 1357 sbuf1_i = sbuf1[pa[tmp2]]; 1358 jmax = sbuf1_i[0]; 1359 ct1 = 2*jmax + 1; 1360 ct2 = 0; 1361 rbuf2_i = rbuf2[tmp2]; 1362 rbuf3_i = rbuf3[tmp2]; 1363 if (!ijonly) rbuf4_i = rbuf4[tmp2]; 1364 for (j=1; j<=jmax; j++) { 1365 is_no = sbuf1_i[2*j-1]; 1366 rmap_i = rmap[is_no]; 1367 if (!allcolumns[is_no]) cmap_i = cmap[is_no]; 1368 subc = (Mat_SeqBAIJ*)submats[is_no]->data; 1369 imat_ilen = subc->ilen; 1370 imat_j = subc->j; 1371 imat_i = subc->i; 1372 if (!ijonly) imat_a = subc->a; 1373 max1 = sbuf1_i[2*j]; 1374 for (k=0; k<max1; k++,ct1++) { /* for each recved block row */ 1375 row = sbuf1_i[ct1]; 1376 1377 if (allrows[is_no]) { 1378 row = sbuf1_i[ct1]; 1379 } else { 1380 #if defined(PETSC_USE_CTABLE) 1381 ierr = PetscTableFind(rmap_i,row+1,&row);CHKERRQ(ierr); 1382 row--; 1383 if (row < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"row not found in table"); 1384 #else 1385 row = rmap_i[row]; 1386 #endif 1387 } 1388 ilen = imat_ilen[row]; 1389 mat_i = imat_i[row]; 1390 if (!ijonly) mat_a = imat_a + mat_i*bs2; 1391 mat_j = imat_j + mat_i; 1392 max2 = rbuf2_i[ct1]; 1393 if (!allcolumns[is_no]) { 1394 for (l=0; l<max2; l++,ct2++) { 1395 #if defined(PETSC_USE_CTABLE) 1396 ierr = PetscTableFind(cmap_i,rbuf3_i[ct2]+1,&tcol);CHKERRQ(ierr); 1397 #else 1398 tcol = cmap_i[rbuf3_i[ct2]]; 1399 #endif 1400 if (tcol) { 1401 *mat_j++ = tcol - 1; 1402 if (!ijonly) { 1403 ierr = PetscMemcpy(mat_a,rbuf4_i+ct2*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1404 mat_a += bs2; 1405 } 1406 ilen++; 1407 } 1408 } 1409 } else { /* allcolumns */ 1410 for (l=0; l<max2; l++,ct2++) { 1411 *mat_j++ = rbuf3_i[ct2]; /* same global column index of C */ 1412 if (!ijonly) { 1413 ierr = PetscMemcpy(mat_a,rbuf4_i+ct2*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1414 mat_a += bs2; 1415 } 1416 ilen++; 1417 } 1418 } 1419 imat_ilen[row] = ilen; 1420 } 1421 } 1422 } 1423 1424 if (!iscsorted) { /* sort column indices of the rows */ 1425 MatScalar *work; 1426 1427 ierr = PetscMalloc1(bs2,&work);CHKERRQ(ierr); 1428 for (i=0; i<ismax; i++) { 1429 subc = (Mat_SeqBAIJ*)submats[i]->data; 1430 imat_ilen = subc->ilen; 1431 imat_j = subc->j; 1432 imat_i = subc->i; 1433 if (!ijonly) imat_a = subc->a; 1434 if (allcolumns[i]) continue; 1435 1436 jmax = nrow[i]; 1437 for (j=0; j<jmax; j++) { 1438 mat_i = imat_i[j]; 1439 mat_j = imat_j + mat_i; 1440 ilen = imat_ilen[j]; 1441 if (ijonly) { 1442 ierr = PetscSortInt(ilen,mat_j);CHKERRQ(ierr); 1443 } else { 1444 mat_a = imat_a + mat_i*bs2; 1445 ierr = PetscSortIntWithDataArray(ilen,mat_j,mat_a,bs2*sizeof(MatScalar),work);CHKERRQ(ierr); 1446 } 1447 } 1448 } 1449 ierr = PetscFree(work);CHKERRQ(ierr); 1450 } 1451 1452 if (!ijonly) { 1453 ierr = PetscFree(r_status4);CHKERRQ(ierr); 1454 ierr = PetscFree(r_waits4);CHKERRQ(ierr); 1455 if (nrqr) {ierr = MPI_Waitall(nrqr,s_waits4,s_status4);CHKERRQ(ierr);} 1456 ierr = PetscFree(s_waits4);CHKERRQ(ierr); 1457 ierr = PetscFree(s_status4);CHKERRQ(ierr); 1458 } 1459 1460 /* Restore the indices */ 1461 for (i=0; i<ismax; i++) { 1462 if (!allrows[i]) { 1463 ierr = ISRestoreIndices(isrow[i],irow+i);CHKERRQ(ierr); 1464 } 1465 if (!allcolumns[i]) { 1466 ierr = ISRestoreIndices(iscol[i],icol+i);CHKERRQ(ierr); 1467 } 1468 } 1469 1470 for (i=0; i<ismax; i++) { 1471 ierr = MatAssemblyBegin(submats[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1472 ierr = MatAssemblyEnd(submats[i],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1473 } 1474 1475 /* Destroy allocated memory */ 1476 if (!ismax) { 1477 ierr = PetscFree(pa);CHKERRQ(ierr); 1478 ierr = PetscFree4(sbuf1,ptr,tmp,ctr);CHKERRQ(ierr); 1479 for (i=0; i<nrqr; ++i) { 1480 ierr = PetscFree(sbuf2[i]);CHKERRQ(ierr); 1481 } 1482 for (i=0; i<nrqs; ++i) { 1483 ierr = PetscFree(rbuf3[i]);CHKERRQ(ierr); 1484 } 1485 1486 ierr = PetscFree3(sbuf2,req_size,req_source1);CHKERRQ(ierr); 1487 ierr = PetscFree3(req_source2,rbuf2,rbuf3);CHKERRQ(ierr); 1488 } 1489 ierr = PetscFree5(irow,icol,nrow,ncol,issorted);CHKERRQ(ierr); 1490 ierr = PetscFree6(smats,row2proc,cmap,rmap,allcolumns,allrows);CHKERRQ(ierr); 1491 1492 if (!ijonly) { 1493 ierr = PetscFree(sbuf_aa[0]);CHKERRQ(ierr); 1494 ierr = PetscFree(sbuf_aa);CHKERRQ(ierr); 1495 1496 for (i=0; i<nrqs; ++i) { 1497 ierr = PetscFree(rbuf4[i]);CHKERRQ(ierr); 1498 } 1499 ierr = PetscFree(rbuf4);CHKERRQ(ierr); 1500 } 1501 c->ijonly = PETSC_FALSE; /* set back to the default */ 1502 PetscFunctionReturn(0); 1503 } 1504